From: Peter Munch Date: Tue, 22 Jun 2021 07:40:05 +0000 (+0200) Subject: Introduce MGTwoLevelTransfer::reinit() X-Git-Tag: v9.4.0-rc1~1185^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=refs%2Fpull%2F12484%2Fhead;p=dealii.git Introduce MGTwoLevelTransfer::reinit() --- diff --git a/doc/news/changes/minor/20210623Munch-1 b/doc/news/changes/minor/20210623Munch-1 new file mode 100644 index 0000000000..2acd580142 --- /dev/null +++ b/doc/news/changes/minor/20210623Munch-1 @@ -0,0 +1,5 @@ +New: The new unified function MGTwoLevelTransfer::reinit() selects automatically +if MGTwoLevelTransfer::reinit_geometric_transfer() or +MGTwoLevelTransfer::reinit_polynomial_transfer() is needed. +
+(Peter Munch, 2021/06/23) diff --git a/examples/step-75/step-75.cc b/examples/step-75/step-75.cc index b946d3455a..951d173b01 100644 --- a/examples/step-75/step-75.cc +++ b/examples/step-75/step-75.cc @@ -700,9 +700,8 @@ namespace Step75 fe_index_for_degree[degree] = i; } - unsigned int minlevel = 0; - unsigned int minlevel_p = n_h_levels; - unsigned int maxlevel = n_h_levels + n_p_levels - 1; + unsigned int minlevel = 0; + unsigned int maxlevel = n_h_levels + n_p_levels - 1; dof_handlers.resize(minlevel, maxlevel); operators.resize(minlevel, maxlevel); @@ -802,17 +801,11 @@ namespace Step75 // Set up intergrid operators and collect transfer operators within a single // operator as needed by the Multigrid solver class. - for (unsigned int level = minlevel; level < minlevel_p; ++level) - transfers[level + 1].reinit_geometric_transfer(dof_handlers[level + 1], - dof_handlers[level], - constraints[level + 1], - constraints[level]); - - for (unsigned int level = minlevel_p; level < maxlevel; ++level) - transfers[level + 1].reinit_polynomial_transfer(dof_handlers[level + 1], - dof_handlers[level], - constraints[level + 1], - constraints[level]); + for (unsigned int level = minlevel; level < maxlevel; ++level) + transfers[level + 1].reinit(dof_handlers[level + 1], + dof_handlers[level], + constraints[level + 1], + constraints[level]); MGTransferGlobalCoarsening transfer( transfers, [&](const auto l, auto &vec) { diff --git a/include/deal.II/multigrid/mg_transfer_global_coarsening.h b/include/deal.II/multigrid/mg_transfer_global_coarsening.h index 9187abc515..38790e66e5 100644 --- a/include/deal.II/multigrid/mg_transfer_global_coarsening.h +++ b/include/deal.II/multigrid/mg_transfer_global_coarsening.h @@ -216,6 +216,27 @@ public: const unsigned int mg_level_fine = numbers::invalid_unsigned_int, const unsigned int mg_level_coarse = numbers::invalid_unsigned_int); + /** + * Set up transfer operator between the given DoFHandler objects ( + * @p dof_handler_fine and @p dof_handler_coarse). Depending on the + * underlying Triangulation objects polynomial or geometrical global + * coarsening is performed. + * + * @note While geometric transfer can be only performed on active levels + * (`numbers::invalid_unsigned_int`), polynomial transfers can also be + * performed on coarse-grid levels. + * + * @note The function polynomial_transfer_supported() can be used to + * check if the given polynomial coarsening strategy is supported. + */ + void + reinit(const DoFHandler & dof_handler_fine, + const DoFHandler & dof_handler_coarse, + const AffineConstraints &constraint_fine, + const AffineConstraints &constraint_coarse, + const unsigned int mg_level_fine = numbers::invalid_unsigned_int, + const unsigned int mg_level_coarse = numbers::invalid_unsigned_int); + /** * Check if a fast templated version of the polynomial transfer between * @p fe_degree_fine and @p fe_degree_coarse is available. diff --git a/include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h b/include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h index 83f55a2374..3e1e5e5e4f 100644 --- a/include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h +++ b/include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h @@ -411,6 +411,8 @@ namespace internal const std::function active_fe_index_function; }; + + template class FineDoFHandlerView { @@ -2804,6 +2806,96 @@ MGTwoLevelTransfer>:: +template +void +MGTwoLevelTransfer>::reinit( + const DoFHandler & dof_handler_fine, + const DoFHandler & dof_handler_coarse, + const AffineConstraints &constraint_fine, + const AffineConstraints &constraint_coarse, + const unsigned int mg_level_fine, + const unsigned int mg_level_coarse) +{ + // determine if polynomial transfer can be performed via the following two + // criteria: + // 1) multigrid levels can be only used with polynomial transfer + bool do_polynomial_transfer = + (mg_level_fine != numbers::invalid_unsigned_int) || + (mg_level_coarse != numbers::invalid_unsigned_int); + + // 2) the meshes are identical + if (do_polynomial_transfer == false) + { + const internal::CellIDTranslator cell_id_translator( + dof_handler_fine.get_triangulation().n_global_coarse_cells(), + dof_handler_fine.get_triangulation().n_global_levels()); + + AssertDimension( + dof_handler_fine.get_triangulation().n_global_coarse_cells(), + dof_handler_coarse.get_triangulation().n_global_coarse_cells()); + AssertIndexRange(dof_handler_coarse.get_triangulation().n_global_levels(), + dof_handler_fine.get_triangulation().n_global_levels() + + 1); + + IndexSet is_locally_owned_fine(cell_id_translator.size()); + IndexSet is_locally_owned_coarse(cell_id_translator.size()); + + for (const auto &cell : dof_handler_fine.active_cell_iterators()) + if (!cell->is_artificial() && cell->is_locally_owned()) + is_locally_owned_fine.add_index(cell_id_translator.translate(cell)); + + for (const auto &cell : dof_handler_coarse.active_cell_iterators()) + if (!cell->is_artificial() && cell->is_locally_owned()) + is_locally_owned_coarse.add_index(cell_id_translator.translate(cell)); + + const MPI_Comm communicator = dof_handler_fine.get_communicator(); + + std::vector owning_ranks( + is_locally_owned_coarse.n_elements()); + + Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload + process(is_locally_owned_fine, + is_locally_owned_coarse, + communicator, + owning_ranks, + false); + + Utilities::MPI::ConsensusAlgorithms::Selector< + std::pair, + unsigned int> + consensus_algorithm(process, communicator); + consensus_algorithm.run(); + + bool all_cells_found = true; + + for (unsigned i = 0; i < is_locally_owned_coarse.n_elements(); ++i) + all_cells_found &= (owning_ranks[i] != numbers::invalid_unsigned_int); + + do_polynomial_transfer = + Utilities::MPI::min(static_cast(all_cells_found), + communicator) == 1; + } + + if (do_polynomial_transfer) + internal::MGTwoLevelTransferImplementation::reinit_polynomial_transfer( + dof_handler_fine, + dof_handler_coarse, + constraint_fine, + constraint_coarse, + mg_level_fine, + mg_level_coarse, + *this); + else + internal::MGTwoLevelTransferImplementation::reinit_geometric_transfer( + dof_handler_fine, + dof_handler_coarse, + constraint_fine, + constraint_coarse, + *this); +} + + + template bool MGTwoLevelTransfer>:: diff --git a/tests/multigrid-global-coarsening/mg_transfer_a_01.cc b/tests/multigrid-global-coarsening/mg_transfer_a_01.cc index 621f8e8b8d..5d0a2a4857 100644 --- a/tests/multigrid-global-coarsening/mg_transfer_a_01.cc +++ b/tests/multigrid-global-coarsening/mg_transfer_a_01.cc @@ -107,10 +107,10 @@ do_test(const FiniteElement &fe_fine, const FiniteElement &fe_coarse) // setup transfer operator MGTwoLevelTransfer> transfer; - transfer.reinit_geometric_transfer(dof_handler_fine, - dof_handler_coarse, - constraint_fine, - constraint_coarse); + transfer.reinit(dof_handler_fine, + dof_handler_coarse, + constraint_fine, + constraint_coarse); test_transfer_operator(transfer, dof_handler_fine, dof_handler_coarse); } diff --git a/tests/multigrid-global-coarsening/mg_transfer_a_02.cc b/tests/multigrid-global-coarsening/mg_transfer_a_02.cc index f9ba6d81af..7f624d96ab 100644 --- a/tests/multigrid-global-coarsening/mg_transfer_a_02.cc +++ b/tests/multigrid-global-coarsening/mg_transfer_a_02.cc @@ -94,10 +94,10 @@ do_test(const FiniteElement &fe_fine, const FiniteElement &fe_coarse) // setup transfer operator MGTwoLevelTransfer> transfer; - transfer.reinit_geometric_transfer(dof_handler_fine, - dof_handler_coarse, - constraint_fine, - constraint_coarse); + transfer.reinit(dof_handler_fine, + dof_handler_coarse, + constraint_fine, + constraint_coarse); test_transfer_operator(transfer, dof_handler_fine, dof_handler_coarse); } diff --git a/tests/multigrid-global-coarsening/mg_transfer_a_03.cc b/tests/multigrid-global-coarsening/mg_transfer_a_03.cc index c956650683..230c2ef27f 100644 --- a/tests/multigrid-global-coarsening/mg_transfer_a_03.cc +++ b/tests/multigrid-global-coarsening/mg_transfer_a_03.cc @@ -102,10 +102,10 @@ do_test(const FiniteElement &fe_fine, const FiniteElement &fe_coarse) // setup transfer operator MGTwoLevelTransfer> transfer; - transfer.reinit_geometric_transfer(dof_handler_fine, - dof_handler_coarse, - constraint_fine, - constraint_coarse); + transfer.reinit(dof_handler_fine, + dof_handler_coarse, + constraint_fine, + constraint_coarse); test_transfer_operator(transfer, dof_handler_fine, dof_handler_coarse); } diff --git a/tests/multigrid-global-coarsening/mg_transfer_p_01.cc b/tests/multigrid-global-coarsening/mg_transfer_p_01.cc index 91561c5327..ada0702cfe 100644 --- a/tests/multigrid-global-coarsening/mg_transfer_p_01.cc +++ b/tests/multigrid-global-coarsening/mg_transfer_p_01.cc @@ -92,10 +92,10 @@ do_test(const FiniteElement &fe_fine, const FiniteElement &fe_coarse) // setup transfer operator MGTwoLevelTransfer> transfer; - transfer.reinit_polynomial_transfer(dof_handler_fine, - dof_handler_coarse, - constraint_fine, - constraint_coarse); + transfer.reinit(dof_handler_fine, + dof_handler_coarse, + constraint_fine, + constraint_coarse); test_transfer_operator(transfer, dof_handler_fine, dof_handler_coarse); } diff --git a/tests/multigrid-global-coarsening/mg_transfer_p_02.cc b/tests/multigrid-global-coarsening/mg_transfer_p_02.cc index 0a32b9630c..0acd948bfc 100644 --- a/tests/multigrid-global-coarsening/mg_transfer_p_02.cc +++ b/tests/multigrid-global-coarsening/mg_transfer_p_02.cc @@ -140,10 +140,10 @@ do_test() // setup transfer operator MGTwoLevelTransfer> transfer; - transfer.reinit_polynomial_transfer(dof_handler_fine, - dof_handler_coarse, - constraint_fine, - constraint_coarse); + transfer.reinit(dof_handler_fine, + dof_handler_coarse, + constraint_fine, + constraint_coarse); test_transfer_operator(transfer, dof_handler_fine, dof_handler_coarse); diff --git a/tests/multigrid-global-coarsening/mg_transfer_p_03.cc b/tests/multigrid-global-coarsening/mg_transfer_p_03.cc index aaa1e21d70..c533f890f4 100644 --- a/tests/multigrid-global-coarsening/mg_transfer_p_03.cc +++ b/tests/multigrid-global-coarsening/mg_transfer_p_03.cc @@ -79,10 +79,10 @@ do_test(const FiniteElement &fe_fine, const FiniteElement &fe_coarse) // setup transfer operator MGTwoLevelTransfer> transfer; - transfer.reinit_polynomial_transfer(dof_handler_fine, - dof_handler_coarse, - constraint_fine, - constraint_coarse); + transfer.reinit(dof_handler_fine, + dof_handler_coarse, + constraint_fine, + constraint_coarse); test_transfer_operator(transfer, dof_handler_fine, dof_handler_coarse); } diff --git a/tests/multigrid-global-coarsening/mg_transfer_p_04.cc b/tests/multigrid-global-coarsening/mg_transfer_p_04.cc index 004ee665a9..1dda8054fd 100644 --- a/tests/multigrid-global-coarsening/mg_transfer_p_04.cc +++ b/tests/multigrid-global-coarsening/mg_transfer_p_04.cc @@ -120,10 +120,10 @@ do_test(const FiniteElement &fe_fine, const FiniteElement &fe_coarse) // setup transfer operator MGTwoLevelTransfer> transfer; - transfer.reinit_polynomial_transfer(dof_handler_fine, - dof_handler_coarse, - constraint_fine, - constraint_coarse); + transfer.reinit(dof_handler_fine, + dof_handler_coarse, + constraint_fine, + constraint_coarse); test_transfer_operator(transfer, dof_handler_fine, dof_handler_coarse); } diff --git a/tests/multigrid-global-coarsening/mg_transfer_p_05.cc b/tests/multigrid-global-coarsening/mg_transfer_p_05.cc index f3a7bf8cfa..60acdc1a7c 100644 --- a/tests/multigrid-global-coarsening/mg_transfer_p_05.cc +++ b/tests/multigrid-global-coarsening/mg_transfer_p_05.cc @@ -92,12 +92,12 @@ do_test(const FiniteElement &fe_fine, const FiniteElement &fe_coarse) deallog.push("coarse"); MGTwoLevelTransfer> transfer; - transfer.reinit_polynomial_transfer(dof_handler_fine, - dof_handler_coarse, - constraint_fine, - constraint_coarse, - 0, - 0); + transfer.reinit(dof_handler_fine, + dof_handler_coarse, + constraint_fine, + constraint_coarse, + 0, + 0); test_transfer_operator(transfer, dof_handler_fine, dof_handler_coarse); deallog.pop(); @@ -106,10 +106,10 @@ do_test(const FiniteElement &fe_fine, const FiniteElement &fe_coarse) deallog.push("active"); MGTwoLevelTransfer> transfer; - transfer.reinit_polynomial_transfer(dof_handler_fine, - dof_handler_coarse, - constraint_fine, - constraint_coarse); + transfer.reinit(dof_handler_fine, + dof_handler_coarse, + constraint_fine, + constraint_coarse); test_transfer_operator(transfer, dof_handler_fine, dof_handler_coarse); deallog.pop(); diff --git a/tests/multigrid-global-coarsening/multigrid_a_01.cc b/tests/multigrid-global-coarsening/multigrid_a_01.cc index fbdc709a90..62c04637a3 100644 --- a/tests/multigrid-global-coarsening/multigrid_a_01.cc +++ b/tests/multigrid-global-coarsening/multigrid_a_01.cc @@ -95,10 +95,10 @@ test(const unsigned int n_refinements, // set up transfer operator for (unsigned int l = min_level; l < max_level; ++l) - transfers[l + 1].reinit_geometric_transfer(dof_handlers[l + 1], - dof_handlers[l], - constraints[l + 1], - constraints[l]); + transfers[l + 1].reinit(dof_handlers[l + 1], + dof_handlers[l], + constraints[l + 1], + constraints[l]); MGTransferGlobalCoarsening transfer( transfers, diff --git a/tests/multigrid-global-coarsening/multigrid_p_01.cc b/tests/multigrid-global-coarsening/multigrid_p_01.cc index fc107374fd..044374f6eb 100644 --- a/tests/multigrid-global-coarsening/multigrid_p_01.cc +++ b/tests/multigrid-global-coarsening/multigrid_p_01.cc @@ -123,10 +123,10 @@ test(const unsigned int n_refinements, // set up transfer operator for (unsigned int l = min_level; l < max_level; ++l) - transfers[l + 1].reinit_polynomial_transfer(dof_handlers[l + 1], - dof_handlers[l], - constraints[l + 1], - constraints[l]); + transfers[l + 1].reinit(dof_handlers[l + 1], + dof_handlers[l], + constraints[l + 1], + constraints[l]); MGTransferGlobalCoarsening transfer( transfers, diff --git a/tests/multigrid-global-coarsening/multigrid_p_02.cc b/tests/multigrid-global-coarsening/multigrid_p_02.cc index 2b7a764838..00aa4430ed 100644 --- a/tests/multigrid-global-coarsening/multigrid_p_02.cc +++ b/tests/multigrid-global-coarsening/multigrid_p_02.cc @@ -98,12 +98,12 @@ test(const unsigned int n_refinements, const unsigned int fe_degree_fine) // set up transfer operator for (unsigned int l = min_level; l < max_level; ++l) - transfers[l + 1].reinit_polynomial_transfer(dof_handlers[l + 1], - dof_handlers[l], - constraints[l + 1], - constraints[l], - 0 /*level*/, - 0 /*level*/); + transfers[l + 1].reinit(dof_handlers[l + 1], + dof_handlers[l], + constraints[l + 1], + constraints[l], + 0 /*level*/, + 0 /*level*/); MGTransferGlobalCoarsening transfer( transfers,