renumber_dofs (const std::vector<types::global_dof_index> &new_numbers,
dealii::DoFHandler<dim,spacedim> &dof_handler,
NumberCache & number_cache_current) const;
+ private:
+
};
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
}
/* --------------------- class ParallelSequential ---------------- */
+
template <int dim, int spacedim>
NumberCache
ParallelShared<dim,spacedim>::
//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;
}
//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 ---------------- */
}
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>
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
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
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