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

Remove deprecated DoFTools::extract_boundary_dofs #15578

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
4 changes: 4 additions & 0 deletions doc/news/changes/incompatibilities/20230702DanielArndt
Original file line number Diff line number Diff line change
@@ -0,0 +1,4 @@
Removed: The deprecated overloads for DoFTools::extract_boundary_dofs
returning void have been removed.
<br>
(Daniel Arndt, 2023/07/02)
74 changes: 0 additions & 74 deletions include/deal.II/dofs/dof_tools.h
Original file line number Diff line number Diff line change
Expand Up @@ -1260,67 +1260,6 @@ namespace DoFTools
const BlockMask & component_mask,
std::vector<bool> & selected_dofs);

/**
* Extract all degrees of freedom which are at the boundary and belong to
* specified components of the solution. The function returns its results in
* the last non-default-valued parameter which contains @p true if a degree
* of freedom is at the boundary and belongs to one of the selected
* components, and @p false otherwise.
*
* By specifying the @p boundary_ids variable, you can select which boundary
* indicators the faces have to have on which the degrees of freedom are
* located that shall be extracted. If it is an empty list, then all
* boundary indicators are accepted.
*
* The size of @p component_mask (see
* @ref GlossComponentMask)
* shall equal the number of components in the finite element used by @p
* dof. The size of @p selected_dofs shall equal
* <tt>dof_handler.n_dofs()</tt>. Previous contents of this array are
* overwritten.
*
* Using the usual convention, if a shape function is non-zero in more than
* one component (i.e. it is non-primitive), then the element in the
* component mask is used that corresponds to the first non-zero components.
* Elements in the mask corresponding to later components are ignored.
*
* @deprecated This function will not work for DoFHandler objects that are built
* on a parallel::distributed::Triangulation object. The reasons is that the
* output argument @p selected_dofs has to have a length equal to <i>all</i>
* global degrees of freedom. Consequently, this does not scale to very
* large problems, and this is also why the function is deprecated. If you
* need the functionality of this function for
* parallel triangulations, then you need to use the other
* DoFTools::extract_boundary_dofs() function that returns its information
* via an IndexSet object.
*
* @param[in] dof_handler The object that describes which degrees of freedom
* live on which cell.
* @param[in] component_mask A mask denoting the vector components of the
* finite element that should be considered (see also
* @ref GlossComponentMask).
* @param[out] selected_dofs A vector of booleans that is returned and for
* which
* an element will be @p true if the corresponding index is a
* degree of freedom that is located on the
* boundary (and correspond to the selected vector components and boundary
* indicators, depending on the values of the @p component_mask and @p
* boundary_ids arguments).
* @param[in] boundary_ids If empty, this function extracts the indices of the
* degrees of freedom for all parts of the boundary. If it is a non- empty
* list, then the function only considers boundary faces with the boundary
* indicators listed in this argument.
*
* @see
* @ref GlossBoundaryIndicator "Glossary entry on boundary indicators"
*/
template <int dim, int spacedim>
DEAL_II_DEPRECATED void
extract_boundary_dofs(const DoFHandler<dim, spacedim> & dof_handler,
const ComponentMask & component_mask,
std::vector<bool> & selected_dofs,
const std::set<types::boundary_id> &boundary_ids = {});

/**
* Extract all degrees of freedom which are at the boundary and belong to
* specified components of the solution. The function returns its results in
Expand Down Expand Up @@ -1369,19 +1308,6 @@ namespace DoFTools
const ComponentMask &component_mask = ComponentMask(),
const std::set<types::boundary_id> &boundary_ids = {});

/**
* The same as the previous function, except that it returns its information
* via the third argument.
*
* @deprecated Use the previous function instead.
*/
template <int dim, int spacedim>
DEAL_II_DEPRECATED void
extract_boundary_dofs(const DoFHandler<dim, spacedim> & dof_handler,
const ComponentMask & component_mask,
IndexSet & selected_dofs,
const std::set<types::boundary_id> &boundary_ids = {});

/**
* This function is similar to the extract_boundary_dofs() function but it
* extracts those degrees of freedom whose shape functions are nonzero on at
Expand Down
4 changes: 2 additions & 2 deletions source/dofs/dof_tools.cc
Original file line number Diff line number Diff line change
Expand Up @@ -554,8 +554,8 @@ namespace DoFTools
"This function can not be used with distributed triangulations. "
"See the documentation for more information."));

IndexSet indices;
extract_boundary_dofs(dof_handler, component_mask, indices, boundary_ids);
IndexSet indices =
extract_boundary_dofs(dof_handler, component_mask, boundary_ids);

// clear and reset array by default values
selected_dofs.clear();
Expand Down
14 changes: 5 additions & 9 deletions tests/bits/anna_4.cc
Original file line number Diff line number Diff line change
Expand Up @@ -189,15 +189,12 @@ FindBug<dim>::dirichlet_conditions()
component_mask);


std::vector<bool> fixed_dofs(dof_handler.n_dofs());
const std::set<types::boundary_id> boundary_ids = {0};

// get a list of those boundary DoFs which
// we want to be fixed:
DoFTools::extract_boundary_dofs(dof_handler,
component_mask,
fixed_dofs,
boundary_ids);
const IndexSet fixed_dofs =
DoFTools::extract_boundary_dofs(dof_handler, component_mask, boundary_ids);

// (Primitive) Check if the DoFs
// where adjusted correctly (note
Expand All @@ -207,7 +204,7 @@ FindBug<dim>::dirichlet_conditions()
// component 0 by 0)
for (unsigned int i = 0; i < dof_handler.n_dofs(); ++i)
{
if (fixed_dofs[i] == true)
if (fixed_dofs.is_element(i))
{
AssertThrow(dirichlet_dofs[i] == 0, ExcInternalError());
}
Expand All @@ -226,9 +223,8 @@ FindBug<dim>::dirichlet_conditions()
VectorBoundaryValues<dim>(),
dirichlet_dofs,
component_mask);
for (unsigned int i = 0; i < dof_handler.n_dofs(); ++i)
if (fixed_dofs[i] == true)
deallog << i << ' ' << dirichlet_dofs[i] << std::endl;
for (const types::global_dof_index dof : fixed_dofs)
deallog << dof << ' ' << dirichlet_dofs[dof] << std::endl;
}


Expand Down
25 changes: 12 additions & 13 deletions tests/bits/anna_6.cc
Original file line number Diff line number Diff line change
Expand Up @@ -180,16 +180,15 @@ ImposeBC<dim>::test_extract_boundary_DoFs()
bc_component_select[1] = true;
bc_component_select[2] = false;

std::vector<bool> ned_boundary_dofs(dof_handler.n_dofs());
const std::set<types::boundary_id> boundary_ids = {0};
DoFTools::extract_boundary_dofs(dof_handler,
bc_component_select,
ned_boundary_dofs,
boundary_ids);
const IndexSet ned_boundary_dofs =
DoFTools::extract_boundary_dofs(dof_handler,
bc_component_select,
boundary_ids);


for (unsigned int i = 0; i < dof_handler.n_dofs(); ++i)
if (ned_boundary_dofs[i] == true)
if (ned_boundary_dofs.is_element(i))
boundary_values[i] = 0.;
}

Expand Down Expand Up @@ -218,12 +217,11 @@ ImposeBC<dim>::test_interpolate_BC()
// check
// (the pressure is assumed to be set to 1
// on the boundary)
std::vector<bool> p_boundary_dofs(dof_handler.n_dofs());
const std::set<types::boundary_id> boundary_ids = {0};
DoFTools::extract_boundary_dofs(dof_handler,
bc_component_select,
p_boundary_dofs,
boundary_ids);
const IndexSet p_boundary_dofs =
DoFTools::extract_boundary_dofs(dof_handler,
bc_component_select,
boundary_ids);
for (unsigned int i = 0; i < dof_handler.n_dofs(); ++i)
{
// error: pressure boundary DoF
Expand All @@ -235,8 +233,9 @@ ImposeBC<dim>::test_interpolate_BC()
// nedelec boundary DoF i has
// wrongly been set to some
// value
AssertThrow((p_boundary_dofs[i] && boundary_values[i] == 1.) ||
(!(p_boundary_dofs[i]) && boundary_values[i] != 1.),
AssertThrow((p_boundary_dofs.is_element(i) && boundary_values[i] == 1.) ||
(!(p_boundary_dofs.is_element(i)) &&
boundary_values[i] != 1.),
ExcInternalError());

deallog << boundary_values[i] << ' ';
Expand Down
12 changes: 4 additions & 8 deletions tests/bits/step-11.cc
Original file line number Diff line number Diff line change
Expand Up @@ -104,19 +104,15 @@ LaplaceProblem<dim>::setup_system()
solution.reinit(dof_handler.n_dofs());
system_rhs.reinit(dof_handler.n_dofs());

std::vector<bool> boundary_dofs(dof_handler.n_dofs(), false);
DoFTools::extract_boundary_dofs(dof_handler,
std::vector<bool>(1, true),
boundary_dofs);
const IndexSet boundary_dofs =
DoFTools::extract_boundary_dofs(dof_handler, std::vector<bool>(1, true));

const unsigned int first_boundary_dof =
std::distance(boundary_dofs.begin(),
std::find(boundary_dofs.begin(), boundary_dofs.end(), true));
const unsigned int first_boundary_dof = *boundary_dofs.begin();

mean_value_constraints.clear();
mean_value_constraints.add_line(first_boundary_dof);
for (unsigned int i = first_boundary_dof + 1; i < dof_handler.n_dofs(); ++i)
if (boundary_dofs[i] == true)
if (boundary_dofs.is_element(i))
mean_value_constraints.add_entry(first_boundary_dof, i, -1);
mean_value_constraints.close();

Expand Down
19 changes: 8 additions & 11 deletions tests/dofs/dof_tools_05.cc
Original file line number Diff line number Diff line change
Expand Up @@ -28,13 +28,12 @@ void
check_this(const DoFHandler<dim> &dof_handler)
{
std::vector<bool> component_select(dof_handler.get_fe().n_components(), true);
std::vector<bool> boundary_dofs(dof_handler.n_dofs());
IndexSet boundary_dofs(dof_handler.n_dofs());

// first with all components
{
DoFTools::extract_boundary_dofs(dof_handler,
component_select,
boundary_dofs);
boundary_dofs =
DoFTools::extract_boundary_dofs(dof_handler, component_select);
output_bool_vector(boundary_dofs);
}

Expand All @@ -43,20 +42,18 @@ check_this(const DoFHandler<dim> &dof_handler)
for (unsigned int i = 1; i < component_select.size(); i += 2)
component_select[i] = false;
{
DoFTools::extract_boundary_dofs(dof_handler,
component_select,
boundary_dofs);
boundary_dofs =
DoFTools::extract_boundary_dofs(dof_handler, component_select);
output_bool_vector(boundary_dofs);
}

// third further restrict to
// boundary indicator 0
{
const std::set<types::boundary_id> boundary_ids = {0};
DoFTools::extract_boundary_dofs(dof_handler,
component_select,
boundary_dofs,
boundary_ids);
boundary_dofs = DoFTools::extract_boundary_dofs(dof_handler,
component_select,
boundary_ids);
output_bool_vector(boundary_dofs);
}
}
11 changes: 4 additions & 7 deletions tests/dofs/interpolate_boundary_values_01.cc
Original file line number Diff line number Diff line change
Expand Up @@ -151,15 +151,12 @@ FindBug<dim>::dirichlet_conditions()
component_mask);


std::vector<bool> fixed_dofs(dof_handler.n_dofs());
const std::set<types::boundary_id> boundary_ids = {0};

// get a list of those boundary DoFs which
// we want to be fixed:
DoFTools::extract_boundary_dofs(dof_handler,
component_mask,
fixed_dofs,
boundary_ids);
const IndexSet fixed_dofs =
DoFTools::extract_boundary_dofs(dof_handler, component_mask, boundary_ids);

// (Primitive) Check if the DoFs
// where adjusted correctly (note
Expand All @@ -169,7 +166,7 @@ FindBug<dim>::dirichlet_conditions()
// component 0 by 0)
for (unsigned int i = 0; i < dof_handler.n_dofs(); ++i)
{
if (fixed_dofs[i] == true)
if (fixed_dofs.is_element(i))
{
AssertThrow(dirichlet_dofs[i] == 13, ExcInternalError());
}
Expand All @@ -189,7 +186,7 @@ FindBug<dim>::dirichlet_conditions()
dirichlet_dofs,
component_mask);
for (unsigned int i = 0; i < dof_handler.n_dofs(); ++i)
if (fixed_dofs[i] == true)
if (fixed_dofs.is_element(i))
deallog << i << ' ' << dirichlet_dofs[i] << std::endl;
}

Expand Down
12 changes: 4 additions & 8 deletions tests/hp/step-11.cc
Original file line number Diff line number Diff line change
Expand Up @@ -107,19 +107,15 @@ LaplaceProblem<dim>::setup_system()
solution.reinit(dof_handler.n_dofs());
system_rhs.reinit(dof_handler.n_dofs());

std::vector<bool> boundary_dofs(dof_handler.n_dofs(), false);
DoFTools::extract_boundary_dofs(dof_handler,
std::vector<bool>(1, true),
boundary_dofs);
const IndexSet boundary_dofs =
DoFTools::extract_boundary_dofs(dof_handler, std::vector<bool>(1, true));

const unsigned int first_boundary_dof =
std::distance(boundary_dofs.begin(),
std::find(boundary_dofs.begin(), boundary_dofs.end(), true));
const unsigned int first_boundary_dof = *boundary_dofs.begin();

mean_value_constraints.clear();
mean_value_constraints.add_line(first_boundary_dof);
for (unsigned int i = first_boundary_dof + 1; i < dof_handler.n_dofs(); ++i)
if (boundary_dofs[i] == true)
if (boundary_dofs.is_element(i))
mean_value_constraints.add_entry(first_boundary_dof, i, -1);
mean_value_constraints.close();

Expand Down
12 changes: 4 additions & 8 deletions tests/meshworker/step-11-mesh_loop.cc
Original file line number Diff line number Diff line change
Expand Up @@ -140,14 +140,10 @@ namespace Step11
solution.reinit(dof_handler.n_dofs());
system_rhs.reinit(dof_handler.n_dofs());

std::vector<bool> boundary_dofs(dof_handler.n_dofs(), false);
DoFTools::extract_boundary_dofs(dof_handler,
ComponentMask(),
boundary_dofs);
const IndexSet boundary_dofs =
DoFTools::extract_boundary_dofs(dof_handler, ComponentMask());

const unsigned int first_boundary_dof = std::distance(
boundary_dofs.begin(),
std::find(boundary_dofs.begin(), boundary_dofs.end(), true));
const unsigned int first_boundary_dof = *boundary_dofs.begin();

// Then generate a constraints object with just this one constraint. First
// clear all previous content (which might reside there from the previous
Expand All @@ -160,7 +156,7 @@ namespace Step11
DoFTools::make_hanging_node_constraints(dof_handler, constraints);
constraints.add_line(first_boundary_dof);
for (unsigned int i = first_boundary_dof + 1; i < dof_handler.n_dofs(); ++i)
if (boundary_dofs[i] == true)
if (boundary_dofs.is_element(i))
constraints.add_entry(first_boundary_dof, i, -1);
constraints.close();

Expand Down