Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Move active_fe_index_type to types. #14210

Merged
merged 1 commit into from
Aug 23, 2022
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
5 changes: 5 additions & 0 deletions include/deal.II/base/types.h
Original file line number Diff line number Diff line change
Expand Up @@ -53,6 +53,11 @@ namespace types
*/
#define DEAL_II_VERTEX_INDEX_MPI_TYPE MPI_UINT64_T

/**
* The type in which we store the active and future FE indices.
*/
using fe_index = unsigned short int;

/**
* The type used to denote the global index of degrees of freedom. This
* type is then also used for querying the global *number* of degrees
Expand Down
17 changes: 9 additions & 8 deletions include/deal.II/dofs/dof_handler.h
Original file line number Diff line number Diff line change
Expand Up @@ -25,6 +25,7 @@
#include <deal.II/base/index_set.h>
#include <deal.II/base/iterator_range.h>
#include <deal.II/base/smartpointer.h>
#include <deal.II/base/types.h>

#include <deal.II/dofs/block_info.h>
#include <deal.II/dofs/dof_accessor.h>
Expand Down Expand Up @@ -523,8 +524,10 @@ class DoFHandler : public Subscriptor

/**
* The type in which we store the active FE index.
*
* @deprecated Use types::fe_index instead.
*/
using active_fe_index_type = unsigned short int;
using active_fe_index_type DEAL_II_DEPRECATED = types::fe_index;

/**
* The type in which we store the offsets in the CRS data structures.
Expand All @@ -535,8 +538,8 @@ class DoFHandler : public Subscriptor
* Invalid active FE index which will be used as a default value to determine
* whether a future FE index has been set or not.
*/
static const active_fe_index_type invalid_active_fe_index =
static_cast<active_fe_index_type>(-1);
static const types::fe_index invalid_active_fe_index =
static_cast<types::fe_index>(-1);
Comment on lines +541 to +542
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Why not move this to numbers?

Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I will do this in a follow-up PR. I'll make small changes, that way we can easily identify if some commit breaks something.


/**
* Standard constructor, not initializing any data. After constructing an
Expand Down Expand Up @@ -1489,7 +1492,7 @@ class DoFHandler : public Subscriptor
* of the appropriate position of a cell in the vectors is done via
* hp_object_fe_ptr (CRS scheme).
*/
mutable std::array<std::vector<active_fe_index_type>, dim + 1>
mutable std::array<std::vector<types::fe_index>, dim + 1>
hp_object_fe_indices;

/**
Expand All @@ -1501,15 +1504,13 @@ class DoFHandler : public Subscriptor
* Active FE index of an active cell (identified by level and level index).
* This vector is only used in hp-mode.
*/
mutable std::vector<std::vector<active_fe_index_type>>
hp_cell_active_fe_indices;
mutable std::vector<std::vector<types::fe_index>> hp_cell_active_fe_indices;

/**
* Future FE index of an active cell (identified by level and level index).
* This vector is only used in hp-mode.
*/
mutable std::vector<std::vector<active_fe_index_type>>
hp_cell_future_fe_indices;
mutable std::vector<std::vector<types::fe_index>> hp_cell_future_fe_indices;

/**
* An array to store the indices for level degrees of freedom located at
Expand Down
47 changes: 18 additions & 29 deletions source/dofs/dof_handler.cc
Original file line number Diff line number Diff line change
Expand Up @@ -47,8 +47,7 @@ const unsigned int DoFHandler<dim, spacedim>::default_fe_index;


template <int dim, int spacedim>
const typename DoFHandler<dim, spacedim>::active_fe_index_type
DoFHandler<dim, spacedim>::invalid_active_fe_index;
const types::fe_index DoFHandler<dim, spacedim>::invalid_active_fe_index;


namespace internal
Expand Down Expand Up @@ -1322,9 +1321,6 @@ namespace internal
dof_handler.hp_capability_enabled == true,
(typename DoFHandler<dim, spacedim>::ExcOnlyAvailableWithHP()));

using active_fe_index_type =
typename dealii::DoFHandler<dim, spacedim>::active_fe_index_type;

if (const dealii::parallel::shared::Triangulation<dim, spacedim> *tr =
dynamic_cast<
const dealii::parallel::shared::Triangulation<dim, spacedim>
Expand All @@ -1340,7 +1336,7 @@ namespace internal
// on the other cells to zero. then we add all of these vectors
// up, and because every vector entry has exactly one processor
// that owns it, the sum is correct
std::vector<active_fe_index_type> active_fe_indices(
std::vector<types::fe_index> active_fe_indices(
tr->n_active_cells(), 0u);
for (const auto &cell : dof_handler.active_cell_iterators())
if (cell->is_locally_owned())
Expand Down Expand Up @@ -1376,17 +1372,15 @@ namespace internal
// to have functions that can pack and unpack the data we want to
// transport -- namely, the single unsigned int active_fe_index
// objects
auto pack =
[](const typename dealii::DoFHandler<dim, spacedim>::
active_cell_iterator &cell) -> active_fe_index_type {
auto pack = [](const typename dealii::DoFHandler<dim, spacedim>::
active_cell_iterator &cell) -> types::fe_index {
return cell->active_fe_index();
};

auto unpack =
[&dof_handler](
const typename dealii::DoFHandler<dim, spacedim>::
active_cell_iterator & cell,
const active_fe_index_type active_fe_index) -> void {
auto unpack = [&dof_handler](
const typename dealii::DoFHandler<dim, spacedim>::
active_cell_iterator &cell,
const types::fe_index active_fe_index) -> void {
// we would like to say
// cell->set_active_fe_index(active_fe_index);
// but this is not allowed on cells that are not
Expand All @@ -1397,7 +1391,7 @@ namespace internal
};

GridTools::exchange_cell_data_to_ghosts<
active_fe_index_type,
types::fe_index,
dealii::DoFHandler<dim, spacedim>>(dof_handler, pack, unpack);
}
else
Expand Down Expand Up @@ -1434,15 +1428,12 @@ namespace internal
dof_handler.hp_capability_enabled == true,
(typename DoFHandler<dim, spacedim>::ExcOnlyAvailableWithHP()));

using active_fe_index_type =
typename dealii::DoFHandler<dim, spacedim>::active_fe_index_type;

if (const dealii::parallel::shared::Triangulation<dim, spacedim> *tr =
dynamic_cast<
const dealii::parallel::shared::Triangulation<dim, spacedim>
*>(&dof_handler.get_triangulation()))
{
std::vector<active_fe_index_type> future_fe_indices(
std::vector<types::fe_index> future_fe_indices(
tr->n_active_cells(), 0u);
for (const auto &cell : dof_handler.active_cell_iterators() |
IteratorFilters::LocallyOwnedCell())
Expand All @@ -1467,26 +1458,24 @@ namespace internal
DistributedTriangulationBase<dim, spacedim> *>(
&dof_handler.get_triangulation()))
{
auto pack =
[&dof_handler](
const typename dealii::DoFHandler<dim, spacedim>::
active_cell_iterator &cell) -> active_fe_index_type {
auto pack = [&dof_handler](
const typename dealii::DoFHandler<dim, spacedim>::
active_cell_iterator &cell) -> types::fe_index {
return dof_handler
.hp_cell_future_fe_indices[cell->level()][cell->index()];
};

auto unpack =
[&dof_handler](
const typename dealii::DoFHandler<dim, spacedim>::
active_cell_iterator & cell,
const active_fe_index_type future_fe_index) -> void {
auto unpack = [&dof_handler](
const typename dealii::DoFHandler<dim, spacedim>::
active_cell_iterator &cell,
const types::fe_index future_fe_index) -> void {
dof_handler
.hp_cell_future_fe_indices[cell->level()][cell->index()] =
future_fe_index;
};

GridTools::exchange_cell_data_to_ghosts<
active_fe_index_type,
types::fe_index,
dealii::DoFHandler<dim, spacedim>>(dof_handler, pack, unpack);
}
else
Expand Down
3 changes: 1 addition & 2 deletions source/hp/refinement.cc
Original file line number Diff line number Diff line change
Expand Up @@ -811,8 +811,7 @@ namespace hp

// there can be as many levels in the hierarchy as active FE indices are
// possible
using level_type =
typename DoFHandler<dim, spacedim>::active_fe_index_type;
using level_type = types::fe_index;
const auto invalid_level = static_cast<level_type>(-1);

// map from FE index to level in hierarchy
Expand Down