Skip to content

Commit

Permalink
Generalize DoFHandlerPolicy of p:d:t
Browse files Browse the repository at this point in the history
  • Loading branch information
peterrum committed Aug 16, 2019
1 parent 6731a01 commit 660234c
Show file tree
Hide file tree
Showing 6 changed files with 130 additions and 100 deletions.
4 changes: 4 additions & 0 deletions doc/news/changes/minor/20190814PeterMunch-1
Original file line number Diff line number Diff line change
@@ -0,0 +1,4 @@
New: Generalize DoFHandlerPolicy of p:d:t such that it uses
the new definition of CellId.
<br>
(Peter Munch, 2019/08/14)
18 changes: 18 additions & 0 deletions include/deal.II/distributed/tria.h
Original file line number Diff line number Diff line change
Expand Up @@ -1217,6 +1217,16 @@ namespace parallel
std::vector<bool>
mark_locally_active_vertices_on_level(const int level) const;

virtual unsigned int
translate_coarse_cell_id_to_coarse_cell_index(
const unsigned int coarse_cell_id) const override;

virtual unsigned int
translate_coarse_cell_index_to_coarse_cell_id(
const unsigned int coarse_cell_index) const override;



template <typename>
friend class dealii::internal::DoFHandlerImplementation::Policy::
ParallelDistributed;
Expand Down Expand Up @@ -1384,6 +1394,14 @@ namespace parallel
*/
virtual std::vector<bool>
mark_locally_active_vertices_on_level(const unsigned int level) const;

virtual unsigned int
translate_coarse_cell_id_to_coarse_cell_index(
const unsigned int coarse_cell_id) const override;

virtual unsigned int
translate_coarse_cell_index_to_coarse_cell_id(
const unsigned int coarse_cell_index) const override;
};
} // namespace distributed
} // namespace parallel
Expand Down
41 changes: 41 additions & 0 deletions source/distributed/tria.cc
Original file line number Diff line number Diff line change
Expand Up @@ -4551,6 +4551,26 @@ namespace parallel



template <int dim, int spacedim>
unsigned int
Triangulation<dim, spacedim>::translate_coarse_cell_id_to_coarse_cell_index(
const unsigned int coarse_cell_id) const
{
return p4est_tree_to_coarse_cell_permutation[coarse_cell_id];
}



template <int dim, int spacedim>
unsigned int
Triangulation<dim, spacedim>::translate_coarse_cell_index_to_coarse_cell_id(
const unsigned int coarse_cell_index) const
{
return coarse_cell_to_p4est_tree_permutation[coarse_cell_index];
}



template <int dim, int spacedim>
void
Triangulation<dim, spacedim>::add_periodicity(
Expand Down Expand Up @@ -5047,6 +5067,27 @@ namespace parallel



template <int spacedim>
unsigned int
Triangulation<1, spacedim>::translate_coarse_cell_id_to_coarse_cell_index(
const unsigned int) const
{
Assert(false, ExcNotImplemented());
return 0;
}



template <int spacedim>
unsigned int
Triangulation<1, spacedim>::translate_coarse_cell_index_to_coarse_cell_id(
const unsigned int) const
{
Assert(false, ExcNotImplemented());
return 0;
}


template <int spacedim>
void
Triangulation<1, spacedim>::load(const std::string &, const bool)
Expand Down
77 changes: 23 additions & 54 deletions source/dofs/dof_handler_policy.cc
Original file line number Diff line number Diff line change
Expand Up @@ -4538,11 +4538,7 @@ namespace internal
communicate_mg_ghost_cells(
const typename parallel::distributed::Triangulation<dim, spacedim>
& tria,
DoFHandlerType &dof_handler,
const std::vector<dealii::types::global_dof_index>
&coarse_cell_to_p4est_tree_permutation,
const std::vector<dealii::types::global_dof_index>
&p4est_tree_to_coarse_cell_permutation)
DoFHandlerType &dof_handler)
{
// build list of cells to request for each neighbor
std::set<dealii::types::subdomain_id> level_ghost_owners =
Expand All @@ -4565,7 +4561,7 @@ namespace internal

find_marked_mg_ghost_cells_recursively<dim, spacedim>(
tria,
coarse_cell_to_p4est_tree_permutation[cell->index()],
cell->id().get_coarse_cell_id(),
cell,
p4est_coarse_cell,
neighbor_cell_list);
Expand Down Expand Up @@ -4634,11 +4630,14 @@ namespace internal
c < cell_data_transfer_buffer.tree_indices.size();
++c)
{
const auto temp =
CellId(cell_data_transfer_buffer.tree_indices[c], 0, NULL)
.to_cell(tria);

typename DoFHandlerType::level_cell_iterator cell(
&dof_handler.get_triangulation(),
0,
p4est_tree_to_coarse_cell_permutation
[cell_data_transfer_buffer.tree_indices[c]],
temp->index(),
&dof_handler);

typename dealii::internal::p4est::types<dim>::quadrant
Expand Down Expand Up @@ -4703,12 +4702,12 @@ namespace internal
c < cell_data_transfer_buffer.tree_indices.size();
++c, dofs += 1 + dofs[0])
{
const auto temp =
CellId(cell_data_transfer_buffer.tree_indices[c], 0, NULL)
.to_cell(tria);

typename DoFHandlerType::level_cell_iterator cell(
&tria,
0,
p4est_tree_to_coarse_cell_permutation
[cell_data_transfer_buffer.tree_indices[c]],
&dof_handler);
&tria, 0, temp->index(), &dof_handler);

typename dealii::internal::p4est::types<dim>::quadrant
p4est_coarse_cell;
Expand Down Expand Up @@ -4750,9 +4749,7 @@ namespace internal
void
communicate_mg_ghost_cells(
const typename parallel::distributed::Triangulation<1, spacedim> &,
DoFHandler<1, spacedim> &,
const std::vector<dealii::types::global_dof_index> &,
const std::vector<dealii::types::global_dof_index> &)
DoFHandler<1, spacedim> &)
{
Assert(false, ExcNotImplemented());
}
Expand All @@ -4763,9 +4760,7 @@ namespace internal
void
communicate_mg_ghost_cells(
const typename parallel::distributed::Triangulation<1, spacedim> &,
hp::DoFHandler<1, spacedim> &,
const std::vector<dealii::types::global_dof_index> &,
const std::vector<dealii::types::global_dof_index> &)
hp::DoFHandler<1, spacedim> &)
{
Assert(false, ExcNotImplemented());
}
Expand Down Expand Up @@ -4794,9 +4789,7 @@ namespace internal
void
communicate_dof_indices_on_marked_cells(
const DoFHandler<1, spacedim> &,
const std::map<unsigned int, std::set<dealii::types::subdomain_id>> &,
const std::vector<dealii::types::global_dof_index> &,
const std::vector<dealii::types::global_dof_index> &)
const std::map<unsigned int, std::set<dealii::types::subdomain_id>> &)
{
Assert(false, ExcNotImplemented());
}
Expand All @@ -4807,9 +4800,7 @@ namespace internal
void
communicate_dof_indices_on_marked_cells(
const hp::DoFHandler<1, spacedim> &,
const std::map<unsigned int, std::set<dealii::types::subdomain_id>> &,
const std::vector<dealii::types::global_dof_index> &,
const std::vector<dealii::types::global_dof_index> &)
const std::map<unsigned int, std::set<dealii::types::subdomain_id>> &)
{
Assert(false, ExcNotImplemented());
}
Expand All @@ -4820,9 +4811,7 @@ namespace internal
void
communicate_dof_indices_on_marked_cells(
const DoFHandlerType &dof_handler,
const std::map<unsigned int, std::set<dealii::types::subdomain_id>> &,
const std::vector<dealii::types::global_dof_index> &,
const std::vector<dealii::types::global_dof_index> &)
const std::map<unsigned int, std::set<dealii::types::subdomain_id>> &)
{
# ifndef DEAL_II_WITH_MPI
(void)vertices_with_ghost_neighbors;
Expand Down Expand Up @@ -5171,10 +5160,7 @@ namespace internal
// as explained in the 'distributed' paper, this has to be
// done twice
communicate_dof_indices_on_marked_cells(
*dof_handler,
vertices_with_ghost_neighbors,
triangulation->coarse_cell_to_p4est_tree_permutation,
triangulation->p4est_tree_to_coarse_cell_permutation);
*dof_handler, vertices_with_ghost_neighbors);

// in case of hp::DoFHandlers, we may have received valid
// indices of degrees of freedom that are dominated by a fe
Expand All @@ -5189,10 +5175,7 @@ namespace internal
// may still have invalid ones. thus, exchange
// one more time.
communicate_dof_indices_on_marked_cells(
*dof_handler,
vertices_with_ghost_neighbors,
triangulation->coarse_cell_to_p4est_tree_permutation,
triangulation->p4est_tree_to_coarse_cell_permutation);
*dof_handler, vertices_with_ghost_neighbors);

// at this point, we must have taken care of the data transfer
// on all cells we had previously marked. verify this
Expand Down Expand Up @@ -5399,11 +5382,7 @@ namespace internal
// Phase 1. Request all marked cells from corresponding owners. If we
// managed to get every DoF, remove the user_flag, otherwise we
// will request them again in the step below.
communicate_mg_ghost_cells(
*triangulation,
*dof_handler,
triangulation->coarse_cell_to_p4est_tree_permutation,
triangulation->p4est_tree_to_coarse_cell_permutation);
communicate_mg_ghost_cells(*triangulation, *dof_handler);

// have a barrier so that sends from above and below this
// place are not mixed up.
Expand All @@ -5423,11 +5402,7 @@ namespace internal

// Phase 2, only request the cells that were not completed
// in Phase 1.
communicate_mg_ghost_cells(
*triangulation,
*dof_handler,
triangulation->coarse_cell_to_p4est_tree_permutation,
triangulation->p4est_tree_to_coarse_cell_permutation);
communicate_mg_ghost_cells(*triangulation, *dof_handler);

# ifdef DEBUG
// make sure we have removed all flags:
Expand Down Expand Up @@ -5708,10 +5683,7 @@ namespace internal
// as explained in the 'distributed' paper, this has to be
// done twice
communicate_dof_indices_on_marked_cells(
*dof_handler,
vertices_with_ghost_neighbors,
triangulation->coarse_cell_to_p4est_tree_permutation,
triangulation->p4est_tree_to_coarse_cell_permutation);
*dof_handler, vertices_with_ghost_neighbors);

// in case of hp::DoFHandlers, we may have received valid
// indices of degrees of freedom that are dominated by a fe
Expand All @@ -5722,10 +5694,7 @@ namespace internal
*dof_handler);

communicate_dof_indices_on_marked_cells(
*dof_handler,
vertices_with_ghost_neighbors,
triangulation->coarse_cell_to_p4est_tree_permutation,
triangulation->p4est_tree_to_coarse_cell_permutation);
*dof_handler, vertices_with_ghost_neighbors);

triangulation->load_user_flags(user_flags);
}
Expand Down

0 comments on commit 660234c

Please sign in to comment.