From: Peter Munch Date: Sat, 26 Jun 2021 17:12:45 +0000 (+0200) Subject: Move CellIDTranslator into a new file X-Git-Tag: v9.4.0-rc1~1190^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=bb4832b03f74cbe7e0d4a926428499c6bf09fddd;p=dealii.git Move CellIDTranslator into a new file --- diff --git a/include/deal.II/grid/cell_id_translator.h b/include/deal.II/grid/cell_id_translator.h new file mode 100644 index 0000000000..ae5b6d8803 --- /dev/null +++ b/include/deal.II/grid/cell_id_translator.h @@ -0,0 +1,254 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2021 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. +// +// --------------------------------------------------------------------- + +#ifndef dealii_cell_id_translator_h +#define dealii_cell_id_translator_h + +#include +#include + +#include + +DEAL_II_NAMESPACE_OPEN + +namespace internal +{ + /** + * A class that helps to translate between CellIDs and globally uniquely + * indices. The resulting index space might be non-contiguous in the case + * of locally refined meshes. + * + * The following code snippet how to set up an IndexSet object for active + * cells, using this class: + * @code + * // set up translator + * CellIDTranslator translator(tria.n_global_coarse_cells(), + * tria.n_global_levels()); + * + * // initialize the index set with the correct size + * IndexSet index_set(translator.size()); + * + * // convert IDs to indices and add them to the index set + * for(const auto & cell : tria.active_cell_iterators()) + * index_set.add_index(translator.translate(cell)); + * @endcode + * + * @note The definition of the indices are not relevant. + */ + template + class CellIDTranslator + { + public: + /** + * Constructor taking the number of global coarse cells and number of + * global levels of a triangulation. + */ + CellIDTranslator(const types::global_cell_index n_coarse_cells, + const types::global_cell_index n_global_levels); + + /** + * Return maximum number of cells, i.e., in the case the triangulation + * would be refined globally. + */ + types::global_cell_index + size() const; + + /** + * Convert a @p cell of type CellAccessor or DoFCellAccessor to an index. + */ + template + types::global_cell_index + translate(const T &cell) const; + + /** + * Convert the @p i-th @p cell of type CellAccessor or DoFCellAccessor to + * an index. + */ + template + types::global_cell_index + translate(const T &cell, const types::global_cell_index i) const; + + /** + * Convert an index to a CellId. + */ + CellId + to_cell_id(const types::global_cell_index id) const; + + private: + /** + * Convert a binary representation to a index. + */ + static types::global_cell_index + convert_cell_id_binary_type_to_level_coarse_cell_id( + const typename CellId::binary_type &binary_representation); + + /** + * Number of global coarse cells. + */ + const types::global_cell_index n_coarse_cells; + + /** + * Number of global levels. + */ + const types::global_cell_index n_global_levels; + + /** + * Offset of each level in the index space. + */ + std::vector tree_sizes; + }; + + + + template + CellIDTranslator::CellIDTranslator( + const types::global_cell_index n_coarse_cells, + const types::global_cell_index n_global_levels) + : n_coarse_cells(n_coarse_cells) + , n_global_levels(n_global_levels) + { + tree_sizes.push_back(0); + for (unsigned int i = 0; i < n_global_levels; ++i) + tree_sizes.push_back( + tree_sizes.back() + + Utilities::pow(GeometryInfo::max_children_per_cell, i) * + n_coarse_cells); + } + + + + template + types::global_cell_index + CellIDTranslator::size() const + { + return n_coarse_cells * + (Utilities::pow(GeometryInfo::max_children_per_cell, + n_global_levels) - + 1); + } + + + + template + template + types::global_cell_index + CellIDTranslator::translate(const T &cell) const + { + types::global_cell_index id = 0; + + id += convert_cell_id_binary_type_to_level_coarse_cell_id( + cell->id().template to_binary()); + + id += tree_sizes[cell->level()]; + + return id; + } + + + + template + template + types::global_cell_index + CellIDTranslator::translate(const T & cell, + const types::global_cell_index i) const + { + return (translate(cell) - tree_sizes[cell->level()]) * + GeometryInfo::max_children_per_cell + + i + tree_sizes[cell->level() + 1]; + } + + + + template + CellId + CellIDTranslator::to_cell_id(const types::global_cell_index id) const + { + std::vector child_indices; + + types::global_cell_index id_temp = id; + + types::global_cell_index level = 0; + + for (; level < n_global_levels; ++level) + if (id < tree_sizes[level]) + break; + level -= 1; + + id_temp -= tree_sizes[level]; + + for (types::global_cell_index l = 0; l < level; ++l) + { + child_indices.push_back(id_temp % + GeometryInfo::max_children_per_cell); + id_temp /= GeometryInfo::max_children_per_cell; + } + + std::reverse(child_indices.begin(), child_indices.end()); + + return {id_temp, child_indices}; // TODO + } + + + + template + types::global_cell_index + CellIDTranslator::convert_cell_id_binary_type_to_level_coarse_cell_id( + const typename CellId::binary_type &binary_representation) + { + // exploiting the structure of CellId::binary_type + // see also the documentation of CellId + + // actual coarse-grid id + const unsigned int coarse_cell_id = binary_representation[0]; + const unsigned int n_child_indices = binary_representation[1] >> 2; + + const unsigned int children_per_value = + sizeof(CellId::binary_type::value_type) * 8 / dim; + unsigned int child_level = 0; + unsigned int binary_entry = 2; + + // path to the get to the cell + std::vector cell_indices; + while (child_level < n_child_indices) + { + Assert(binary_entry < binary_representation.size(), ExcInternalError()); + + for (unsigned int j = 0; j < children_per_value; ++j) + { + unsigned int cell_index = + (((binary_representation[binary_entry] >> (j * dim))) & + (GeometryInfo::max_children_per_cell - 1)); + cell_indices.push_back(cell_index); + ++child_level; + if (child_level == n_child_indices) + break; + } + ++binary_entry; + } + + // compute new coarse-grid id: c_{i+1} = c_{i}*2^dim + q; + types::global_cell_index level_coarse_cell_id = coarse_cell_id; + for (auto i : cell_indices) + level_coarse_cell_id = + level_coarse_cell_id * GeometryInfo::max_children_per_cell + i; + + return level_coarse_cell_id; + } + +} // namespace internal + +DEAL_II_NAMESPACE_CLOSE + +#endif 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 58ef8e9192..f4bec5b795 100644 --- a/include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h +++ b/include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h @@ -32,6 +32,7 @@ #include +#include #include #include @@ -326,52 +327,6 @@ namespace internal return cell_local_chilren_indices; } - template - types::coarse_cell_id - convert_cell_id_binary_type_to_level_coarse_cell_id( - const typename CellId::binary_type &binary_representation) - { - // exploiting the structure of CellId::binary_type - // see also the documentation of CellId - - // actual coarse-grid id - const unsigned int coarse_cell_id = binary_representation[0]; - const unsigned int n_child_indices = binary_representation[1] >> 2; - - const unsigned int children_per_value = - sizeof(CellId::binary_type::value_type) * 8 / dim; - unsigned int child_level = 0; - unsigned int binary_entry = 2; - - // path to the get to the cell - std::vector cell_indices; - while (child_level < n_child_indices) - { - Assert(binary_entry < binary_representation.size(), - ExcInternalError()); - - for (unsigned int j = 0; j < children_per_value; ++j) - { - unsigned int cell_index = - (((binary_representation[binary_entry] >> (j * dim))) & - (GeometryInfo::max_children_per_cell - 1)); - cell_indices.push_back(cell_index); - ++child_level; - if (child_level == n_child_indices) - break; - } - ++binary_entry; - } - - // compute new coarse-grid id: c_{i+1} = c_{i}*2^dim + q; - types::coarse_cell_id level_coarse_cell_id = coarse_cell_id; - for (auto i : cell_indices) - level_coarse_cell_id = - level_coarse_cell_id * GeometryInfo::max_children_per_cell + i; - - return level_coarse_cell_id; - } - template std::unique_ptr> create_1D_fe(const FiniteElement &fe) @@ -456,89 +411,6 @@ namespace internal const std::function active_fe_index_function; }; - template - class CellIDTranslator - { - public: - CellIDTranslator(const types::global_cell_index n_coarse_cells, - const types::global_cell_index n_global_levels) - : n_coarse_cells(n_coarse_cells) - , n_global_levels(n_global_levels) - { - tree_sizes.push_back(0); - for (unsigned int i = 0; i < n_global_levels; ++i) - tree_sizes.push_back( - tree_sizes.back() + - Utilities::pow(GeometryInfo::max_children_per_cell, i) * - n_coarse_cells); - } - - types::global_cell_index - size() const - { - return n_coarse_cells * - (Utilities::pow(GeometryInfo::max_children_per_cell, - n_global_levels) - - 1); - } - - template - types::global_cell_index - translate(const T &cell) const - { - types::global_cell_index id = 0; - - id += convert_cell_id_binary_type_to_level_coarse_cell_id( - cell->id().template to_binary()); - - id += tree_sizes[cell->level()]; - - return id; - } - - template - types::global_cell_index - translate(const T &cell, const types::global_cell_index i) const - { - return (translate(cell) - tree_sizes[cell->level()]) * - GeometryInfo::max_children_per_cell + - i + tree_sizes[cell->level() + 1]; - } - - CellId - to_cell_id(const types::global_cell_index id) const - { - std::vector child_indices; - - types::global_cell_index id_temp = id; - - types::global_cell_index level = 0; - - for (; level < n_global_levels; ++level) - if (id < tree_sizes[level]) - break; - level -= 1; - - id_temp -= tree_sizes[level]; - - for (types::global_cell_index l = 0; l < level; ++l) - { - child_indices.push_back(id_temp % - GeometryInfo::max_children_per_cell); - id_temp /= GeometryInfo::max_children_per_cell; - } - - std::reverse(child_indices.begin(), child_indices.end()); - - return {id_temp, child_indices}; // TODO - } - - private: - const types::global_cell_index n_coarse_cells; - const types::global_cell_index n_global_levels; - std::vector tree_sizes; - }; - template class FineDoFHandlerView { diff --git a/source/distributed/repartitioning_policy_tools.cc b/source/distributed/repartitioning_policy_tools.cc index a7e35036f4..b814dcf591 100644 --- a/source/distributed/repartitioning_policy_tools.cc +++ b/source/distributed/repartitioning_policy_tools.cc @@ -18,8 +18,9 @@ #include #include +#include -#include +#include DEAL_II_NAMESPACE_OPEN