From: Timo Heister Date: Fri, 23 Oct 2015 15:22:34 +0000 (-0400) Subject: optimize and simplify computing relevant_dofs X-Git-Tag: v8.4.0-rc2~277^2~4 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=691ccdbdcb47546c3dc9ada35acc1347feda8d71;p=dealii.git optimize and simplify computing relevant_dofs simplify logic in extract_locally_relevant_dofs and extract_locally_relevant_mg_dofs. A benchmark shows that the new approach is faster. --- diff --git a/source/dofs/dof_tools.cc b/source/dofs/dof_tools.cc index 55ea3f8664..7ae323db0d 100644 --- a/source/dofs/dof_tools.cc +++ b/source/dofs/dof_tools.cc @@ -946,11 +946,14 @@ namespace DoFTools // collect all the locally owned dofs dof_set = dof_handler.locally_owned_dofs(); - // add the DoF on the adjacent ghost cells to the IndexSet, cache them - // in a set. need to check each dof manually because we can't be sure - // that the dof range of locally_owned_dofs is really contiguous. + // now add the DoF on the adjacent ghost cells to the IndexSet + + // Note: It is not worth it to cache intermediate data in a set because + // add_indices is more efficient. I also benchmarked skipping the + // locally_owned_dofs by doing locally_owned_dofs().is_element() but + // that is also slower unless 70%+ of the DoFs are locally owned and they + // are contiguous. - Timo Heister std::vector dof_indices; - std::set global_dof_indices; typename DH::active_cell_iterator cell = dof_handler.begin_active(), endc = dof_handler.end(); @@ -959,16 +962,9 @@ namespace DoFTools { dof_indices.resize(cell->get_fe().dofs_per_cell); cell->get_dof_indices(dof_indices); - - for (std::vector::iterator it=dof_indices.begin(); - it!=dof_indices.end(); - ++it) - if (!dof_set.is_element(*it)) - global_dof_indices.insert(*it); + dof_set.add_indices(dof_indices.begin(), dof_indices.end()); } - dof_set.add_indices(global_dof_indices.begin(), global_dof_indices.end()); - dof_set.compress(); } @@ -982,11 +978,14 @@ namespace DoFTools // collect all the locally owned dofs dof_set = dof_handler.locally_owned_mg_dofs(level); - // add the DoF on the adjacent ghost cells to the IndexSet, cache them - // in a set. need to check each dof manually because we can't be sure - // that the dof range of locally_owned_dofs is really contiguous. + // add the DoF on the adjacent ghost cells to the IndexSet + + // Note: It is not worth it to cache intermediate data in a set because + // add_indices is more efficient. I also benchmarked skipping the + // locally_owned_dofs by doing locally_owned_dofs().is_element() but + // that is also slower unless 70%+ of the DoFs are locally owned and they + // are contiguous. - Timo Heister std::vector dof_indices; - std::set global_dof_indices; typename DH::cell_iterator cell = dof_handler.begin(level), endc = dof_handler.end(level); @@ -1001,16 +1000,9 @@ namespace DoFTools dof_indices.resize(cell->get_fe().dofs_per_cell); cell->get_mg_dof_indices(dof_indices); - - for (std::vector::iterator it=dof_indices.begin(); - it!=dof_indices.end(); - ++it) - if (!dof_set.is_element(*it)) - global_dof_indices.insert(*it); + dof_set.add_indices(dof_indices.begin(), dof_indices.end()); } - dof_set.add_indices(global_dof_indices.begin(), global_dof_indices.end()); - dof_set.compress(); }