From: marcfehling Date: Fri, 23 Aug 2019 00:35:40 +0000 (+0200) Subject: hp::DoFHandler: Provide own CoarseningStrategy while transferring active_fe_indices... X-Git-Tag: v9.2.0-rc1~1177^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=refs%2Fpull%2F8637%2Fhead;p=dealii.git hp::DoFHandler: Provide own CoarseningStrategy while transferring active_fe_indices on p::d::Triangulations. --- diff --git a/source/hp/dof_handler.cc b/source/hp/dof_handler.cc index c4b84a7b47..db2a329660 100644 --- a/source/hp/dof_handler.cc +++ b/source/hp/dof_handler.cc @@ -37,8 +37,6 @@ #include #include -#include - #include #include @@ -1228,6 +1226,40 @@ namespace internal cell->set_active_fe_index(coarsen.second); } } + + + /** + * Coarsening strategy for the CellDataTransfer object responsible for + * tranferring the active_fe_index of each cell on + * parallel::distributed::Triangulation objects that have been refined. + * + * A finite element index needs to be determined for the (not yet + * active) parent cell from its (still active) children. Out of the set + * of elements previously assigned to the former children, we choose the + * one dominated by all children for the parent cell. + */ + template + static unsigned int + determine_fe_from_children( + const std::vector & children_fe_indices, + dealii::hp::FECollection &fe_collection) + { + Assert(!children_fe_indices.empty(), ExcInternalError()); + + // convert vector to set + const std::set children_fe_indices_set( + children_fe_indices.begin(), children_fe_indices.end()); + + const unsigned int dominated_fe_index = + fe_collection.find_dominated_fe_extended(children_fe_indices_set, + /*codim=*/0); + + Assert(dominated_fe_index != numbers::invalid_unsigned_int, + typename dealii::hp::FECollection< + dim>::ExcNoDominatedFiniteElementAmongstChildren()); + + return dominated_fe_index; + } }; } // namespace DoFHandlerImplementation } // namespace hp @@ -2041,45 +2073,21 @@ namespace hp active_fe_index_transfer = std_cxx14::make_unique(); - // First, do what we would do in the sequential case. - dealii::internal::hp::DoFHandlerImplementation::Implementation:: - collect_fe_indices_on_cells_to_be_refined(*this); - // If we work on a p::d::Triangulation, we have to transfer all // active_fe_indices since ownership of cells may change. We will // use our p::d::CellDataTransfer member to achieve this. Further, // we prepare the values in such a way that they will correspond to // the active_fe_indices on the new mesh. - // Gather all current active_fe_indices. - get_active_fe_indices(active_fe_index_transfer->active_fe_indices); - - // Overwrite active_fe_indices of cells that change by either - // h/p refinement/coarsening. - for (const auto &pair : - active_fe_index_transfer->persisting_cells_fe_index) - active_fe_index_transfer - ->active_fe_indices[pair.first->active_cell_index()] = pair.second; - - for (const auto &pair : - active_fe_index_transfer->refined_cells_fe_index) - active_fe_index_transfer - ->active_fe_indices[pair.first->active_cell_index()] = pair.second; - - for (const auto &pair : - active_fe_index_transfer->coarsened_cells_fe_index) - for (unsigned int child_index = 0; - child_index < pair.first->n_children(); - ++child_index) - { - // Make sure that all children belong to the same subdomain. - Assert(pair.first->child(child_index)->is_locally_owned(), - ExcInternalError()); + // Gather all current future_fe_indices. + active_fe_index_transfer->active_fe_indices.resize( + get_triangulation().n_active_cells(), numbers::invalid_unsigned_int); - active_fe_index_transfer - ->active_fe_indices[pair.first->child(child_index) - ->active_cell_index()] = pair.second; - } + for (const auto &cell : active_cell_iterators()) + if (cell->is_locally_owned()) + active_fe_index_transfer + ->active_fe_indices[cell->active_cell_index()] = + cell->future_fe_index(); // Create transfer object and attach to it. const auto *distributed_tria = dynamic_cast< @@ -2091,7 +2099,10 @@ namespace hp CellDataTransfer>>( *distributed_tria, /*transfer_variable_size_data=*/false, - &CoarseningStrategies::check_equality); + std::bind(&dealii::internal::hp::DoFHandlerImplementation:: + Implementation::determine_fe_from_children, + std::placeholders::_1, + std::ref(fe_collection))); active_fe_index_transfer->cell_data_transfer ->prepare_for_coarsening_and_refinement( @@ -2198,7 +2209,10 @@ namespace hp CellDataTransfer>>( *distributed_tria, /*transfer_variable_size_data=*/false, - &CoarseningStrategies::check_equality); + std::bind(&dealii::internal::hp::DoFHandlerImplementation:: + Implementation::determine_fe_from_children, + std::placeholders::_1, + std::ref(fe_collection))); // If we work on a p::d::Triangulation, we have to transfer all // active fe indices since ownership of cells may change. @@ -2271,7 +2285,10 @@ namespace hp CellDataTransfer>>( *distributed_tria, /*transfer_variable_size_data=*/false, - &CoarseningStrategies::check_equality); + std::bind(&dealii::internal::hp::DoFHandlerImplementation:: + Implementation::determine_fe_from_children, + std::placeholders::_1, + std::ref(fe_collection))); // Unpack active_fe_indices. active_fe_index_transfer->active_fe_indices.resize( diff --git a/tests/mpi/solution_transfer_05.cc b/tests/mpi/solution_transfer_05.cc new file mode 100644 index 0000000000..cdc4a8aaff --- /dev/null +++ b/tests/mpi/solution_transfer_05.cc @@ -0,0 +1,127 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 1998 - 2018 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, redistribute +// it, and/or modify it under the terms of the GNU Lesser General +// Public License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// The full text of the license can be found in the file LICENSE.md at +// the top level directory of deal.II. +// +// --------------------------------------------------------------------- + + + +#include + +#include +#include + +#include +#include + +#include + +#include + +#include +#include + +#include + +#include + +#include +#include + +#include "../tests.h" + + +template +void +transfer(std::ostream &out) +{ + parallel::distributed::Triangulation tria(MPI_COMM_WORLD); + GridGenerator::hyper_cube(tria); + tria.refine_global(2); + + const unsigned int max_degree = 6 - dim; + hp::FECollection fe_dgq; + for (unsigned int deg = 1; deg <= max_degree; ++deg) + fe_dgq.push_back(FE_Q(deg)); + + hp::DoFHandler dgq_dof_handler(tria); + + // randomly assign FE orders + for (const auto &cell : dgq_dof_handler.active_cell_iterators()) + if (cell->is_locally_owned()) + cell->set_active_fe_index(Testing::rand() % max_degree); + dgq_dof_handler.distribute_dofs(fe_dgq); + + IndexSet dgq_locally_owned_dofs = dgq_dof_handler.locally_owned_dofs(); + IndexSet dgq_locally_relevant_dofs; + dealii::DoFTools::extract_locally_relevant_dofs(dgq_dof_handler, + dgq_locally_relevant_dofs); + IndexSet dgq_ghost_dofs = dgq_locally_relevant_dofs; + dgq_ghost_dofs.subtract_set(dgq_locally_owned_dofs); + + LinearAlgebra::distributed::Vector dgq_solution; + dgq_solution.reinit(dgq_locally_owned_dofs, dgq_ghost_dofs, MPI_COMM_WORLD); + + VectorTools::interpolate(dgq_dof_handler, ZeroFunction(), dgq_solution); + + parallel::distributed::SolutionTransfer< + dim, + LinearAlgebra::distributed::Vector, + hp::DoFHandler> + dgq_soltrans(dgq_dof_handler); + + LinearAlgebra::distributed::Vector dgq_old_solution = dgq_solution; + dgq_old_solution.update_ghost_values(); + { + unsigned int counter = 0; + for (auto cell = tria.begin_active(); cell != tria.end(); ++cell, ++counter) + if (cell->is_locally_owned()) + { + if (counter > ((dim == 2) ? 4 : 8)) + cell->set_coarsen_flag(); + else + cell->set_refine_flag(); + } + } + + tria.prepare_coarsening_and_refinement(); + dgq_soltrans.prepare_for_coarsening_and_refinement(dgq_old_solution); + tria.execute_coarsening_and_refinement(); + + dgq_dof_handler.distribute_dofs(fe_dgq); + + dgq_locally_owned_dofs = dgq_dof_handler.locally_owned_dofs(); + dealii::DoFTools::extract_locally_relevant_dofs(dgq_dof_handler, + dgq_locally_relevant_dofs); + dgq_ghost_dofs = dgq_locally_relevant_dofs; + dgq_ghost_dofs.subtract_set(dgq_locally_owned_dofs); + + dgq_solution.reinit(dgq_locally_owned_dofs, dgq_ghost_dofs, MPI_COMM_WORLD); + dgq_soltrans.interpolate(dgq_solution); + + deallog << "OK" << std::endl; +} + + +int +main(int argc, char *argv[]) +{ + Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1); + MPILogInitAll log; + + deallog.push("2d"); + transfer<2>(deallog.get_file_stream()); + deallog.pop(); + deallog.push("3d"); + transfer<3>(deallog.get_file_stream()); + deallog.pop(); +} diff --git a/tests/mpi/solution_transfer_05.with_p4est=true.with_trilinos=true.mpirun=1.output b/tests/mpi/solution_transfer_05.with_p4est=true.with_trilinos=true.mpirun=1.output new file mode 100644 index 0000000000..995f6e446b --- /dev/null +++ b/tests/mpi/solution_transfer_05.with_p4est=true.with_trilinos=true.mpirun=1.output @@ -0,0 +1,3 @@ + +DEAL:0:2d::OK +DEAL:0:3d::OK diff --git a/tests/mpi/solution_transfer_05.with_p4est=true.with_trilinos=true.mpirun=2.output b/tests/mpi/solution_transfer_05.with_p4est=true.with_trilinos=true.mpirun=2.output new file mode 100644 index 0000000000..4aa906268b --- /dev/null +++ b/tests/mpi/solution_transfer_05.with_p4est=true.with_trilinos=true.mpirun=2.output @@ -0,0 +1,7 @@ + +DEAL:0:2d::OK +DEAL:0:3d::OK + +DEAL:1:2d::OK +DEAL:1:3d::OK + diff --git a/tests/mpi/solution_transfer_05.with_p4est=true.with_trilinos=true.mpirun=4.output b/tests/mpi/solution_transfer_05.with_p4est=true.with_trilinos=true.mpirun=4.output new file mode 100644 index 0000000000..c0650f6365 --- /dev/null +++ b/tests/mpi/solution_transfer_05.with_p4est=true.with_trilinos=true.mpirun=4.output @@ -0,0 +1,15 @@ + +DEAL:0:2d::OK +DEAL:0:3d::OK + +DEAL:1:2d::OK +DEAL:1:3d::OK + + +DEAL:2:2d::OK +DEAL:2:3d::OK + + +DEAL:3:2d::OK +DEAL:3:3d::OK +