--- /dev/null
+New: The new unified function MGTwoLevelTransfer::reinit() selects automatically
+if MGTwoLevelTransfer::reinit_geometric_transfer() or
+MGTwoLevelTransfer::reinit_polynomial_transfer() is needed.
+<br>
+(Peter Munch, 2021/06/23)
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);
// 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<dim, VectorType> transfer(
transfers, [&](const auto l, auto &vec) {
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<dim> & dof_handler_fine,
+ const DoFHandler<dim> & dof_handler_coarse,
+ const AffineConstraints<Number> &constraint_fine,
+ const AffineConstraints<Number> &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.
const std::function<unsigned int()> active_fe_index_function;
};
+
+
template <int dim>
class FineDoFHandlerView
{
+template <int dim, typename Number>
+void
+MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::reinit(
+ const DoFHandler<dim> & dof_handler_fine,
+ const DoFHandler<dim> & dof_handler_coarse,
+ const AffineConstraints<Number> &constraint_fine,
+ const AffineConstraints<Number> &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<dim> 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<unsigned int> 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<types::global_cell_index, types::global_cell_index>,
+ 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<unsigned int>(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 <int dim, typename Number>
bool
MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
// setup transfer operator
MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>> 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);
}
// setup transfer operator
MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>> 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);
}
// setup transfer operator
MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>> 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);
}
// setup transfer operator
MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>> 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);
}
// setup transfer operator
MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>
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);
// setup transfer operator
MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>> 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);
}
// setup transfer operator
MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>> 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.push("coarse");
MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>
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();
deallog.push("active");
MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>
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();
// 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<dim, VectorType> transfer(
transfers,
// 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<dim, VectorType> transfer(
transfers,
// 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<dim, VectorType> transfer(
transfers,