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

Change get_active_fe_indices() to return the result. #14207

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 doc/news/changes/incompatibilities/20220820Fehling
Original file line number Diff line number Diff line change
@@ -0,0 +1,5 @@
Changed: The function DoFHandler::get_active_fe_indices() now returns
the result vector. The previous version that writes into the argument
as well as DoFTools::get_active_fe_indices() have been deprecated.
<br>
(Marc Fehling, 08/20/2022)
7 changes: 1 addition & 6 deletions include/deal.II/dofs/dof_accessor.templates.h
Original file line number Diff line number Diff line change
Expand Up @@ -2333,15 +2333,10 @@ namespace internal
accessor.dof_handler->hp_cell_future_fe_indices.size(),
ExcMessage("DoFHandler not initialized"));

// TODO
using active_fe_index_type = unsigned short int;
static const active_fe_index_type invalid_active_fe_index =
static_cast<active_fe_index_type>(-1);

accessor.dof_handler
->hp_cell_future_fe_indices[accessor.level()]
[accessor.present_index] =
invalid_active_fe_index;
DoFHandler<dim, spacedim>::invalid_active_fe_index;
}
};
} // namespace DoFCellAccessorImplementation
Expand Down
16 changes: 13 additions & 3 deletions include/deal.II/dofs/dof_handler.h
Original file line number Diff line number Diff line change
Expand Up @@ -573,16 +573,26 @@ class DoFHandler : public Subscriptor

/**
* Go through the triangulation and set the active FE indices of all
* active cells to the values given in @p active_fe_indices.
* locally owned cells to the values given in @p active_fe_indices.
*/
void
set_active_fe_indices(const std::vector<unsigned int> &active_fe_indices);

/**
* Go through the triangulation and return a vector of active FE indices of
* all locally relevant cells.
*/
std::vector<unsigned int>
get_active_fe_indices() const;

/**
* Go through the triangulation and store the active FE indices of all
* active cells to the vector @p active_fe_indices. This vector is
* resized, if necessary.
* locally relevant cells to the vector @p active_fe_indices. This vector
* is resized, if necessary.
*
* @deprecated Use the function that returns the result vector.
*/
DEAL_II_DEPRECATED_EARLY
void
get_active_fe_indices(std::vector<unsigned int> &active_fe_indices) const;

Expand Down
4 changes: 3 additions & 1 deletion include/deal.II/dofs/dof_tools.h
Original file line number Diff line number Diff line change
Expand Up @@ -2133,9 +2133,11 @@ namespace DoFTools
* For DoFHandler objects without hp-capabilities given as first argument, the
* returned vector will consist of only zeros, indicating that all cells use
* the same finite element. In hp-mode, the values may be different, though.
*
* @deprecated Use DoFHandler::get_active_fe_indices() instead.
*/
template <int dim, int spacedim>
void
DEAL_II_DEPRECATED_EARLY void
get_active_fe_indices(const DoFHandler<dim, spacedim> &dof_handler,
std::vector<unsigned int> & active_fe_indices);

Expand Down
24 changes: 18 additions & 6 deletions source/dofs/dof_handler.cc
Original file line number Diff line number Diff line change
@@ -1,6 +1,6 @@
// ---------------------------------------------------------------------
//
// Copyright (C) 1998 - 2021 by the deal.II authors
Copy link
Member

Choose a reason for hiding this comment

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

This happens automatically before release ;)

// Copyright (C) 1998 - 2022 by the deal.II authors
//
// This file is part of the deal.II library.
//
Expand Down Expand Up @@ -2645,18 +2645,30 @@ DoFHandler<dim, spacedim>::set_active_fe_indices(


template <int dim, int spacedim>
void
DoFHandler<dim, spacedim>::get_active_fe_indices(
std::vector<unsigned int> &active_fe_indices) const
std::vector<unsigned int>
DoFHandler<dim, spacedim>::get_active_fe_indices() const
{
active_fe_indices.resize(this->get_triangulation().n_active_cells());
std::vector<unsigned int> active_fe_indices(
this->get_triangulation().n_active_cells(), invalid_active_fe_index);

// we could try to extract the values directly, since they are
// stored as protected data of this object, but for simplicity we
// use the cell-wise access.
for (const auto &cell : this->active_cell_iterators())
if (!cell->is_artificial())
active_fe_indices[cell->active_cell_index()] = cell->active_fe_index();

return active_fe_indices;
}



template <int dim, int spacedim>
void
DoFHandler<dim, spacedim>::get_active_fe_indices(
std::vector<unsigned int> &active_fe_indices) const
{
active_fe_indices = get_active_fe_indices();
}


Expand Down Expand Up @@ -3026,7 +3038,7 @@ DoFHandler<dim, spacedim>::prepare_for_serialization_of_active_fe_indices()
// active FE indices since ownership of cells may change.

// Gather all current active FE indices
get_active_fe_indices(active_fe_index_transfer->active_fe_indices);
active_fe_index_transfer->active_fe_indices = get_active_fe_indices();

// Attach to transfer object
active_fe_index_transfer->cell_data_transfer->prepare_for_serialization(
Expand Down
3 changes: 1 addition & 2 deletions source/dofs/dof_tools.cc
Original file line number Diff line number Diff line change
Expand Up @@ -1373,8 +1373,7 @@ namespace DoFTools
AssertDimension(active_fe_indices.size(),
dof_handler.get_triangulation().n_active_cells());

for (const auto &cell : dof_handler.active_cell_iterators())
active_fe_indices[cell->active_cell_index()] = cell->active_fe_index();
active_fe_indices = dof_handler.get_active_fe_indices();
}

template <int dim, int spacedim>
Expand Down
13 changes: 5 additions & 8 deletions tests/hp/get_active_fe_indices.cc
Original file line number Diff line number Diff line change
@@ -1,6 +1,6 @@
// ---------------------------------------------------------------------
//
// Copyright (C) 2005 - 2020 by the deal.II authors
// Copyright (C) 2005 - 2022 by the deal.II authors
//
// This file is part of the deal.II library.
//
Expand All @@ -16,7 +16,7 @@


// distribute different finite elements randomly across the domain, then use
// DoFTools::get_active_fe_indices()
// DoFHandler::get_active_fe_indices()


#include <deal.II/dofs/dof_accessor.h>
Expand Down Expand Up @@ -54,16 +54,13 @@ test()

DoFHandler<dim> dof_handler(tria);

for (typename DoFHandler<dim>::active_cell_iterator cell =
dof_handler.begin_active();
cell != dof_handler.end();
++cell)
for (const auto &cell : dof_handler.active_cell_iterators())
cell->set_active_fe_index(Testing::rand() % fe_collection.size());

dof_handler.distribute_dofs(fe_collection);

std::vector<unsigned int> active_fe_indices(tria.n_active_cells());
DoFTools::get_active_fe_indices(dof_handler, active_fe_indices);
std::vector<unsigned int> active_fe_indices =
dof_handler.get_active_fe_indices();
for (unsigned int i = 0; i < tria.n_active_cells(); ++i)
deallog << active_fe_indices[i] << std::endl;
}
Expand Down