]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
reworked DoFTools::locally_owned_dofs_with_subdomain; use updated version in policy
authordenis.davydov <denis.davydov@0785d39b-7218-0410-832d-ea1e28bc413d>
Fri, 2 May 2014 09:40:49 +0000 (09:40 +0000)
committerdenis.davydov <denis.davydov@0785d39b-7218-0410-832d-ea1e28bc413d>
Fri, 2 May 2014 09:40:49 +0000 (09:40 +0000)
git-svn-id: https://svn.dealii.org/branches/branch_sharedtria@32872 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/source/dofs/dof_handler_policy.cc
deal.II/source/dofs/dof_tools.cc
deal.II/source/dofs/dof_tools.inst.in

index 9456dc491a9f0a8d283266906a0c302415a9c6ef..03a56e22ff81a2dc50822c88d04741661b577d01 100644 (file)
@@ -177,6 +177,8 @@ namespace internal
          renumber_dofs (const std::vector<types::global_dof_index>  &new_numbers,
                         dealii::DoFHandler<dim,spacedim> &dof_handler,
                         NumberCache & number_cache_current) const;
+      private:
+
       };
 
 
index 12f13864c7fa661bb29078539b176e50fa7bcc7b..b08b1f12f62be29f2977e73750effc7e1e4d8377 100644 (file)
@@ -1373,15 +1373,13 @@ namespace DoFTools
                                  IndexSet &dof_set);
 
   /**
-   * Extract the set of global DoF indices distributed on
-   * cells with the given subdomain id. This function should be used for
+   * Extract the set of global DoF indices distributed on cells. This function should be used for
    * DoFHandler objects built on dealii::Triangulation or parallel::shared::Triangulation
    * classes. Internally the funciton is based on DoFTools::get_subdomain_association.
    */
   template <class DH>
-  IndexSet
-  locally_owned_dofs_with_subdomain (const DH                  &dof_handler,
-                                            const types::subdomain_id  subdomain   );
+  std::vector<IndexSet>
+  locally_owned_dofs_with_subdomain (const DH   &dof_handler);
 
   /**
    * For each DoF, return in the output array to which subdomain (as
index c5d0aef91d740e195736dc809f296728e2f9b212..0e03d7342b185d789ed9ae6b8839d3c27d44a1dc 100644 (file)
@@ -966,6 +966,7 @@ namespace internal
       }
 
       /* --------------------- class ParallelSequential ---------------- */
+
       template <int dim, int spacedim>
       NumberCache
       ParallelShared<dim,spacedim>::
@@ -976,12 +977,11 @@ namespace internal
          //update current number cache in DoFHandler
          number_cache_current = number_cache;
          DoFRenumbering::subdomain_wise (dof_handler);
-         const parallel::shared::Triangulation<dim,spacedim> *shared_tr = dynamic_cast<const parallel::shared::Triangulation<dim,spacedim>*>(&dof_handler.get_tria());
-         const unsigned int n_mpi_processes =  Utilities::MPI::n_mpi_processes( shared_tr->get_communicator () );
-         number_cache.n_locally_owned_dofs_per_processor.resize (n_mpi_processes);
-         for (unsigned int i=0; i<n_mpi_processes; ++i)
-                 number_cache.n_locally_owned_dofs_per_processor[i]  = DoFTools::count_dofs_with_subdomain_association (dof_handler, i);
-         number_cache.locally_owned_dofs = dealii::DoFTools::locally_owned_dofs_with_subdomain(dof_handler,dof_handler.get_tria().locally_owned_subdomain() );
+         number_cache.locally_owned_dofs_per_processor = DoFTools::locally_owned_dofs_with_subdomain (dof_handler);
+         number_cache.locally_owned_dofs = number_cache.locally_owned_dofs_per_processor[dof_handler.get_tria().locally_owned_subdomain()];
+         number_cache.n_locally_owned_dofs_per_processor.resize (number_cache.locally_owned_dofs_per_processor.size());
+         for (unsigned int i = 0; i < number_cache.n_locally_owned_dofs_per_processor.size(); i++)
+                 number_cache.n_locally_owned_dofs_per_processor[i] = number_cache.locally_owned_dofs_per_processor[i].n_elements();
          return number_cache;
       }
 
@@ -1006,8 +1006,9 @@ namespace internal
           //update current number cache
           number_cache_current = number_cache;
          //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;
+                 number_cache.locally_owned_dofs_per_processor = DoFTools::locally_owned_dofs_with_subdomain (dof_handler);
+                 number_cache.locally_owned_dofs = number_cache.locally_owned_dofs_per_processor[dof_handler.get_tria().locally_owned_subdomain()];
+                 return number_cache;
       }
 
       /* --------------------- class ParallelDistributed ---------------- */
index d45f4b48d12a90bafb544e6cbbf1a839699b77f2..5f01c410d2234582afac52a1a3e5634b751d553a 100644 (file)
@@ -1076,42 +1076,65 @@ namespace DoFTools
   }
 
   template <class DH>
-  IndexSet
-  locally_owned_dofs_with_subdomain (const DH                  &dof_handler,
-                                            const types::subdomain_id  subdomain   )
+  std::vector<IndexSet>
+  locally_owned_dofs_with_subdomain (const DH  &dof_handler)
   {
-         IndexSet index_set (dof_handler.n_dofs());
+         //the following is a random process (flip of a coin), thus should be called once only.
          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;
+
+         const unsigned int num_subdomains = 1 + (*max_element(subdomain_association.begin(), subdomain_association.end()));
+
+         std::vector<IndexSet> index_sets (num_subdomains,IndexSet(dof_handler.n_dofs()));
+         std::vector<types::global_dof_index> i_min(num_subdomains,0),
+                                                      i_max_plus_1(num_subdomains,0);
 
          //loop over subdomain_association
          //and populate IndexSet when
          //a given subdomain is found
-         bool found_start = false;
+         std::vector<bool> found_start(num_subdomains,false);
          for (types::global_dof_index ind = 0; ind < subdomain_association.size(); ind++)
          {
-                 if (!found_start && (subdomain_association[ind]== subdomain) )
+                 for (unsigned int i = 0; i < num_subdomains; i++)
                  {
-                         found_start = true;
-                         i_min = ind;
+                         if (!found_start[i] && (subdomain_association[ind]== i) )
+                         {
+                                 found_start[i] = true;
+                                 i_min[i] = ind;
+                                 break;
+                         }
+                         if (found_start[i] && (subdomain_association[ind] != i) )
+                         {
+                                 found_start[i] = false;
+                                 i_max_plus_1[i] = ind;
+                                 index_sets[i].add_range(i_min[i],i_max_plus_1[i]);
+                                 break;
+                         }
                  }
-                 if (found_start && (subdomain_association[ind] != subdomain) )
+         }
+         //check that only 1 found_start[i] is true:
+         {
+                 bool found = false;
+                 for (unsigned int i = 0; i < num_subdomains; i++)
                  {
-                         found_start = false;
-                         i_max_plus_1 = ind;
-                         index_set.add_range(i_min,i_max_plus_1);
+                         if (found_start[i])
+                         {
+                                 Assert(!found, ExcInternalError())
+                                 found = true;
+                         }
                  }
          }
-
          // 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());
+         for (unsigned int i = 0;i< num_subdomains; i++)
+         {
+                 if (found_start[i]) {
+                         index_sets[i].add_range(i_min[i],subdomain_association.size());
+                 }
+                 index_sets[i].compress();
          }
-         index_set.compress();
-         return index_set;
+         return index_sets;
   }
 
   template <class DH>
index 4e16d7b6d878e8268d36785d8c0512c7fe81c6c0..b0e66feaf85e48f5ea895a9b36e05c1b8110ba57 100644 (file)
@@ -271,15 +271,13 @@ DoFTools::get_subdomain_association<hp::DoFHandler<deal_II_dimension> >
  std::vector<types::subdomain_id>           &subdomain_association);
  
 template
-IndexSet
+std::vector<IndexSet>
 DoFTools::locally_owned_dofs_with_subdomain<DoFHandler<deal_II_dimension> >
-(const DoFHandler<deal_II_dimension> &dof_handler,
- const types::subdomain_id            subdomain);
+(const DoFHandler<deal_II_dimension> &dof_handler);
 template
-IndexSet
+std::vector<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);
+(const hp::DoFHandler<deal_II_dimension> &dof_handler);
 
 template
 unsigned int
@@ -364,15 +362,13 @@ DoFTools::get_subdomain_association<hp::DoFHandler<deal_II_dimension,deal_II_dim
  std::vector<types::subdomain_id>           &subdomain_association);
  
 template
-IndexSet
+std::vector<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);
+(const DoFHandler<deal_II_dimension,deal_II_dimension+1> &dof_handler);
 template
-IndexSet
+std::vector<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);
+(const hp::DoFHandler<deal_II_dimension,deal_II_dimension+1> &dof_handler);
 
 template
 void
@@ -403,15 +399,13 @@ DoFTools::get_subdomain_association<hp::DoFHandler<1,3> >
  std::vector<types::subdomain_id>           &subdomain_association);
  
 template
-IndexSet
+std::vector<IndexSet>
 DoFTools::locally_owned_dofs_with_subdomain<DoFHandler<1,3> >
-(const DoFHandler<1,3>     &dof_handler,
- const types::subdomain_id  subdomain);
+(const DoFHandler<1,3>     &dof_handler);
 template
-IndexSet
+std::vector<IndexSet>
 DoFTools::locally_owned_dofs_with_subdomain<hp::DoFHandler<1,3> >
-(const hp::DoFHandler<1,3> &dof_handler,
- const types::subdomain_id  subdomain);
+(const hp::DoFHandler<1,3> &dof_handler);
 
 template
 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.