From: Martin Kronbichler Date: Mon, 1 Nov 2021 18:20:53 +0000 (+0100) Subject: Simplify code by using cell->is_locally_owned() and friends X-Git-Tag: v9.4.0-rc1~868^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=refs%2Fpull%2F12904%2Fhead;p=dealii.git Simplify code by using cell->is_locally_owned() and friends --- diff --git a/include/deal.II/matrix_free/matrix_free.templates.h b/include/deal.II/matrix_free/matrix_free.templates.h index eafa6b1afa..3a7dc80fc1 100644 --- a/include/deal.II/matrix_free/matrix_free.templates.h +++ b/include/deal.II/matrix_free/matrix_free.templates.h @@ -835,14 +835,12 @@ namespace internal template void resolve_cell(const InIterator & cell, - std::vector> &cell_its, - const unsigned int subdomain_id) + std::vector> &cell_its) { if (cell->has_children()) for (unsigned int child = 0; child < cell->n_children(); ++child) - resolve_cell(cell->child(child), cell_its, subdomain_id); - else if (subdomain_id == numbers::invalid_subdomain_id || - cell->subdomain_id() == subdomain_id) + resolve_cell(cell->child(child), cell_its); + else if (cell->is_locally_owned()) { Assert(cell->is_active(), ExcInternalError()); cell_its.emplace_back(cell->level(), cell->index()); @@ -868,32 +866,20 @@ MatrixFree::initialize_dof_handlers( for (unsigned int no = 0; no < dof_handlers.size(); ++no) dof_info[no].vectorization_length = VectorizedArrayType::size(); - const unsigned int n_mpi_procs = task_info.n_procs; - const unsigned int my_pid = task_info.my_pid; - const Triangulation &tria = dof_handlers[0]->get_triangulation(); const unsigned int level = additional_data.mg_level; if (level == numbers::invalid_unsigned_int) { - if (n_mpi_procs == 1) - cell_level_index.reserve(tria.n_active_cells()); - // For serial Triangulations always take all cells - const unsigned int subdomain_id = - (dynamic_cast *>( - &dof_handlers[0]->get_triangulation()) != nullptr) ? - my_pid : - numbers::invalid_subdomain_id; + cell_level_index.reserve(tria.n_active_cells()); // Go through cells on zeroth level and then successively step down into // children. This gives a z-ordering of the cells, which is beneficial // when setting up neighboring relations between cells for thread // parallelization for (const auto &cell : tria.cell_iterators_on_level(0)) - internal::MatrixFreeFunctions::resolve_cell(cell, - cell_level_index, - subdomain_id); + internal::MatrixFreeFunctions::resolve_cell(cell, cell_level_index); - Assert(n_mpi_procs > 1 || + Assert(task_info.n_procs > 1 || cell_level_index.size() == tria.n_active_cells(), ExcInternalError()); } @@ -904,7 +890,7 @@ MatrixFree::initialize_dof_handlers( { cell_level_index.reserve(tria.n_cells(level)); for (const auto &cell : tria.cell_iterators_on_level(level)) - if (cell->level_subdomain_id() == my_pid) + if (cell->is_locally_owned_on_level()) cell_level_index.emplace_back(cell->level(), cell->index()); } } diff --git a/source/multigrid/mg_tools.cc b/source/multigrid/mg_tools.cc index 02f4fb12b1..706f3b24cd 100644 --- a/source/multigrid/mg_tools.cc +++ b/source/multigrid/mg_tools.cc @@ -606,13 +606,8 @@ namespace MGTools const unsigned int dofs_per_cell = dof.get_fe().n_dofs_per_cell(); std::vector dofs_on_this_cell(dofs_per_cell); - typename DoFHandler::cell_iterator cell = dof.begin(level), - endc = dof.end(level); - for (; cell != endc; ++cell) - if (dof.get_triangulation().locally_owned_subdomain() == - numbers::invalid_subdomain_id || - cell->level_subdomain_id() == - dof.get_triangulation().locally_owned_subdomain()) + for (const auto &cell : dof.cell_iterators_on_level(level)) + if (cell->is_locally_owned_on_level()) { cell->get_mg_dof_indices(dofs_on_this_cell); constraints.add_entries_local_to_global(dofs_on_this_cell, @@ -1348,9 +1343,7 @@ namespace MGTools { for (const auto &cell : dof.cell_iterators()) { - if (dof.get_triangulation().locally_owned_subdomain() != - numbers::invalid_subdomain_id && - cell->level_subdomain_id() == numbers::artificial_subdomain_id) + if (cell->is_artificial_on_level()) continue; const FiniteElement &fe = cell->get_fe(); const unsigned int level = cell->level(); @@ -1380,9 +1373,7 @@ namespace MGTools "It's probably worthwhile to select at least one component.")); for (const auto &cell : dof.cell_iterators()) - if (dof.get_triangulation().locally_owned_subdomain() == - numbers::invalid_subdomain_id || - cell->level_subdomain_id() != numbers::artificial_subdomain_id) + if (!cell->is_artificial_on_level()) for (const unsigned int face_no : GeometryInfo::face_indices()) { if (cell->at_boundary(face_no) == false) @@ -1508,9 +1499,7 @@ namespace MGTools { // Do not look at artificial level cells (in a serial computation we // need to ignore the level_subdomain_id() because it is never set). - if (mg_dof_handler.get_triangulation().locally_owned_subdomain() != - numbers::invalid_subdomain_id && - cell->level_subdomain_id() == numbers::artificial_subdomain_id) + if (cell->is_artificial_on_level()) continue; bool has_coarser_neighbor = false; @@ -1529,11 +1518,7 @@ namespace MGTools // only process cell pairs if one or both of them are owned by // me (ignore if running in serial) - if (mg_dof_handler.get_triangulation() - .locally_owned_subdomain() != - numbers::invalid_subdomain_id && - neighbor->level_subdomain_id() == - numbers::artificial_subdomain_id) + if (neighbor->is_artificial_on_level()) continue; // Do refinement face from the coarse side diff --git a/source/multigrid/mg_transfer_prebuilt.cc b/source/multigrid/mg_transfer_prebuilt.cc index e986774fac..5897d5b3c0 100644 --- a/source/multigrid/mg_transfer_prebuilt.cc +++ b/source/multigrid/mg_transfer_prebuilt.cc @@ -215,17 +215,11 @@ MGTransferPrebuilt::build( DoFTools::extract_locally_relevant_level_dofs(dof_handler, level + 1, level_p1_relevant_dofs); - DynamicSparsityPattern dsp(this->sizes[level + 1], + DynamicSparsityPattern dsp(this->sizes[level + 1], this->sizes[level], level_p1_relevant_dofs); - typename DoFHandler::cell_iterator cell, - endc = dof_handler.end(level); - for (cell = dof_handler.begin(level); cell != endc; ++cell) - if (cell->has_children() && - (dof_handler.get_triangulation().locally_owned_subdomain() == - numbers::invalid_subdomain_id || - cell->level_subdomain_id() == - dof_handler.get_triangulation().locally_owned_subdomain())) + for (const auto &cell : dof_handler.cell_iterators_on_level(level)) + if (cell->has_children() && cell->is_locally_owned_on_level()) { cell->get_mg_dof_indices(dof_indices_parent); @@ -301,12 +295,8 @@ MGTransferPrebuilt::build( FullMatrix prolongation; // now actually build the matrices - for (cell = dof_handler.begin(level); cell != endc; ++cell) - if (cell->has_children() && - (dof_handler.get_triangulation().locally_owned_subdomain() == - numbers::invalid_subdomain_id || - cell->level_subdomain_id() == - dof_handler.get_triangulation().locally_owned_subdomain())) + for (const auto &cell : dof_handler.cell_iterators_on_level(level)) + if (cell->has_children() && cell->is_locally_owned_on_level()) { cell->get_mg_dof_indices(dof_indices_parent);