Skip to content

Commit

Permalink
Merge pull request #8302 from mathsen/extract_dofs_parallel
Browse files Browse the repository at this point in the history
Changed extract_dofs so that it works in parallel
  • Loading branch information
masterleinad committed Jun 11, 2019
2 parents f6544a1 + e04d12b commit 30eb8db
Show file tree
Hide file tree
Showing 62 changed files with 1,633 additions and 228 deletions.
3 changes: 3 additions & 0 deletions doc/news/changes/minor/20190607ArndtAnselmann
Original file line number Diff line number Diff line change
@@ -0,0 +1,3 @@
Unified DoFTools::extract_dofs() for DoFHandler and hp::DoFHandler.
<br>
(Daniel Arndt, Mathias Anselmann, 2019/06/07)
39 changes: 12 additions & 27 deletions include/deal.II/dofs/dof_tools.h
Original file line number Diff line number Diff line change
Expand Up @@ -1346,22 +1346,16 @@ namespace DoFTools
* corresponds to a vector component selected by the mask above. The size
* of this array must equal DoFHandler::n_locally_owned_dofs(), which for
* sequential computations of course equals DoFHandler::n_dofs(). The
* previous contents of this array are overwritten.
*/
template <int dim, int spacedim>
void
extract_dofs(const DoFHandler<dim, spacedim> &dof_handler,
const ComponentMask & component_mask,
std::vector<bool> & selected_dofs);

/**
* The same function as above, but for a hp::DoFHandler.
* previous contents of this array are overwritten. Note that the resulting
* vector just holds the locally owned extracted degrees of freedom, which
* first have to be mapped to the global degrees of freedom, to correspond
* with them.
*/
template <int dim, int spacedim>
template <typename DoFHandlerType>
void
extract_dofs(const hp::DoFHandler<dim, spacedim> &dof_handler,
const ComponentMask & component_mask,
std::vector<bool> & selected_dofs);
extract_dofs(const DoFHandlerType &dof_handler,
const ComponentMask & component_mask,
std::vector<bool> & selected_dofs);

/**
* This function is the equivalent to the DoFTools::extract_dofs() functions
Expand All @@ -1387,20 +1381,11 @@ namespace DoFTools
* sequential computations of course equals DoFHandler::n_dofs(). The
* previous contents of this array are overwritten.
*/
template <int dim, int spacedim>
void
extract_dofs(const DoFHandler<dim, spacedim> &dof_handler,
const BlockMask & block_mask,
std::vector<bool> & selected_dofs);

/**
* The same function as above, but for a hp::DoFHandler.
*/
template <int dim, int spacedim>
template <typename DoFHandlerType>
void
extract_dofs(const hp::DoFHandler<dim, spacedim> &dof_handler,
const BlockMask & block_mask,
std::vector<bool> & selected_dofs);
extract_dofs(const DoFHandlerType &dof_handler,
const BlockMask & block_mask,
std::vector<bool> & selected_dofs);

/**
* Do the same thing as the corresponding extract_dofs() function for one
Expand Down
94 changes: 18 additions & 76 deletions source/dofs/dof_tools.cc
Original file line number Diff line number Diff line change
Expand Up @@ -397,16 +397,14 @@ namespace DoFTools



template <int dim, int spacedim>
template <typename DoFHandlerType>
void
extract_dofs(const DoFHandler<dim, spacedim> &dof,
const ComponentMask & component_mask,
std::vector<bool> & selected_dofs)
extract_dofs(const DoFHandlerType &dof,
const ComponentMask & component_mask,
std::vector<bool> & selected_dofs)
{
const FiniteElement<dim, spacedim> &fe = dof.get_fe();
(void)fe;

Assert(component_mask.represents_n_components(fe.n_components()),
Assert(component_mask.represents_n_components(
dof.get_fe_collection().n_components()),
ExcMessage(
"The given component mask is not sized correctly to represent the "
"components of the given finite element."));
Expand All @@ -416,13 +414,15 @@ namespace DoFTools

// two special cases: no component is selected, and all components are
// selected; both rather stupid, but easy to catch
if (component_mask.n_selected_components(n_components(dof)) == 0)
if (component_mask.n_selected_components(
dof.get_fe_collection().n_components()) == 0)
{
std::fill_n(selected_dofs.begin(), dof.n_locally_owned_dofs(), false);
return;
}
else if (component_mask.n_selected_components(n_components(dof)) ==
n_components(dof))
else if (component_mask.n_selected_components(
dof.get_fe_collection().n_components()) ==
dof.get_fe_collection().n_components())
{
std::fill_n(selected_dofs.begin(), dof.n_locally_owned_dofs(), true);
return;
Expand All @@ -433,7 +433,7 @@ namespace DoFTools
std::fill_n(selected_dofs.begin(), dof.n_locally_owned_dofs(), false);

// get the component association of each DoF and then select the ones
// that match the given set of blocks
// that match the given set of components
std::vector<unsigned char> dofs_by_component(dof.n_locally_owned_dofs());
internal::get_component_association(dof, component_mask, dofs_by_component);

Expand All @@ -443,74 +443,16 @@ namespace DoFTools
}


// TODO: Unify the following two functions with the non-hp case

template <int dim, int spacedim>
void
extract_dofs(const hp::DoFHandler<dim, spacedim> &dof,
const ComponentMask & component_mask,
std::vector<bool> & selected_dofs)
{
const FiniteElement<dim, spacedim> &fe = dof.begin_active()->get_fe();
(void)fe;

Assert(component_mask.represents_n_components(fe.n_components()),
ExcMessage(
"The given component mask is not sized correctly to represent the "
"components of the given finite element."));
Assert(selected_dofs.size() == dof.n_dofs(),
ExcDimensionMismatch(selected_dofs.size(), dof.n_dofs()));

// two special cases: no component is selected, and all components are
// selected; both rather stupid, but easy to catch
if (component_mask.n_selected_components(n_components(dof)) == 0)
{
std::fill_n(selected_dofs.begin(), dof.n_dofs(), false);
return;
}
else if (component_mask.n_selected_components(n_components(dof)) ==
n_components(dof))
{
std::fill_n(selected_dofs.begin(), dof.n_dofs(), true);
return;
}


// preset all values by false
std::fill_n(selected_dofs.begin(), dof.n_dofs(), false);

// get the component association of each DoF and then select the ones
// that match the given set of components
std::vector<unsigned char> dofs_by_component(dof.n_dofs());
internal::get_component_association(dof, component_mask, dofs_by_component);

for (types::global_dof_index i = 0; i < dof.n_dofs(); ++i)
if (component_mask[dofs_by_component[i]] == true)
selected_dofs[i] = true;
}



template <int dim, int spacedim>
void
extract_dofs(const DoFHandler<dim, spacedim> &dof,
const BlockMask & block_mask,
std::vector<bool> & selected_dofs)
{
// simply forward to the function that works based on a component mask
extract_dofs(dof, dof.get_fe().component_mask(block_mask), selected_dofs);
}



template <int dim, int spacedim>
template <typename DoFHandlerType>
void
extract_dofs(const hp::DoFHandler<dim, spacedim> &dof,
const BlockMask & block_mask,
std::vector<bool> & selected_dofs)
extract_dofs(const DoFHandlerType &dof,
const BlockMask & block_mask,
std::vector<bool> & selected_dofs)
{
// simply forward to the function that works based on a component mask
extract_dofs(dof, dof.get_fe().component_mask(block_mask), selected_dofs);
extract_dofs<DoFHandlerType>(
dof, dof.get_fe_collection().component_mask(block_mask), selected_dofs);
}


Expand Down
8 changes: 4 additions & 4 deletions source/dofs/dof_tools.inst.in
Original file line number Diff line number Diff line change
Expand Up @@ -689,25 +689,25 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : DIMENSIONS)
&patch);

template void
extract_dofs<deal_II_dimension, deal_II_space_dimension>(
extract_dofs<DoFHandler<deal_II_dimension, deal_II_space_dimension>>(
const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
const ComponentMask &,
std::vector<bool> &);

template void
extract_dofs<deal_II_dimension, deal_II_space_dimension>(
extract_dofs<DoFHandler<deal_II_dimension, deal_II_space_dimension>>(
const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
const BlockMask &,
std::vector<bool> &);

template void
extract_dofs<deal_II_dimension, deal_II_space_dimension>(
extract_dofs<hp::DoFHandler<deal_II_dimension, deal_II_space_dimension>>(
const hp::DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
const ComponentMask &,
std::vector<bool> &);

template void
extract_dofs<deal_II_dimension, deal_II_space_dimension>(
extract_dofs<hp::DoFHandler<deal_II_dimension, deal_II_space_dimension>>(
const hp::DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
const BlockMask &,
std::vector<bool> &);
Expand Down
4 changes: 3 additions & 1 deletion source/dofs/dof_tools_constraints.cc
Original file line number Diff line number Diff line change
Expand Up @@ -3133,7 +3133,9 @@ namespace DoFTools
{
std::vector<bool> mask(coarse_grid.get_fe(0).n_components(), false);
mask[coarse_component] = true;
extract_dofs(coarse_grid, ComponentMask(mask), coarse_dof_is_parameter);
extract_dofs<DoFHandler<dim, spacedim>>(coarse_grid,
ComponentMask(mask),
coarse_dof_is_parameter);
}

// now we know that the weights in each row constitute a constraint. enter
Expand Down
5 changes: 3 additions & 2 deletions tests/dofs/dof_tools_00a.cc
Original file line number Diff line number Diff line change
Expand Up @@ -16,16 +16,17 @@

#include "../tests.h"
#include "dof_tools_common.h"
#include "dof_tools_common_fake_hp.h"

// check
// DoFTools::
// count_dofs_per_component (...);



template <int dim>
template <typename DoFHandlerType>
void
check_this(const DoFHandler<dim> &dof_handler)
check_this(const DoFHandlerType &dof_handler)
{
const unsigned int n_components = dof_handler.get_fe().n_components();

Expand Down
5 changes: 3 additions & 2 deletions tests/dofs/dof_tools_01a.cc
Original file line number Diff line number Diff line change
Expand Up @@ -18,6 +18,7 @@

#include "../tests.h"
#include "dof_tools_common.h"
#include "dof_tools_common_fake_hp.h"

// check
// DoFTools::
Expand All @@ -26,9 +27,9 @@



template <int dim>
template <typename DoFHandlerType>
void
check_this(const DoFHandler<dim> &dof_handler)
check_this(const DoFHandlerType &dof_handler)
{
// create sparsity pattern
SparsityPattern sp(dof_handler.n_dofs(),
Expand Down
5 changes: 3 additions & 2 deletions tests/dofs/dof_tools_01a_constraints_false.cc
Original file line number Diff line number Diff line change
Expand Up @@ -18,6 +18,7 @@

#include "../tests.h"
#include "dof_tools_common.h"
#include "dof_tools_common_fake_hp.h"

// check
// DoFTools::
Expand All @@ -28,9 +29,9 @@



template <int dim>
template <typename DoFHandlerType>
void
check_this(const DoFHandler<dim> &dof_handler)
check_this(const DoFHandlerType &dof_handler)
{
AffineConstraints<double> cm;
DoFTools::make_hanging_node_constraints(dof_handler, cm);
Expand Down
5 changes: 3 additions & 2 deletions tests/dofs/dof_tools_01a_constraints_true.cc
Original file line number Diff line number Diff line change
Expand Up @@ -18,6 +18,7 @@

#include "../tests.h"
#include "dof_tools_common.h"
#include "dof_tools_common_fake_hp.h"

// check
// DoFTools::
Expand All @@ -28,9 +29,9 @@



template <int dim>
template <typename DoFHandlerType>
void
check_this(const DoFHandler<dim> &dof_handler)
check_this(const DoFHandlerType &dof_handler)
{
AffineConstraints<double> cm;
DoFTools::make_hanging_node_constraints(dof_handler, cm);
Expand Down
5 changes: 3 additions & 2 deletions tests/dofs/dof_tools_01a_subdomain.cc
Original file line number Diff line number Diff line change
Expand Up @@ -18,6 +18,7 @@

#include "../tests.h"
#include "dof_tools_common.h"
#include "dof_tools_common_fake_hp.h"

// check
// DoFTools::
Expand All @@ -27,9 +28,9 @@



template <int dim>
template <typename DoFHandlerType>
void
check_this(const DoFHandler<dim> &dof_handler)
check_this(const DoFHandlerType &dof_handler)
{
// create sparsity pattern
SparsityPattern sp(dof_handler.n_dofs(),
Expand Down
5 changes: 3 additions & 2 deletions tests/dofs/dof_tools_01b.cc
Original file line number Diff line number Diff line change
Expand Up @@ -18,6 +18,7 @@

#include "../tests.h"
#include "dof_tools_common.h"
#include "dof_tools_common_fake_hp.h"

// check
// DoFTools::
Expand All @@ -26,9 +27,9 @@



template <int dim>
template <typename DoFHandlerType>
void
check_this(const DoFHandler<dim> &dof_handler)
check_this(const DoFHandlerType &dof_handler)
{
// create sparsity pattern
DynamicSparsityPattern sp(dof_handler.n_dofs());
Expand Down
5 changes: 3 additions & 2 deletions tests/dofs/dof_tools_01c.cc
Original file line number Diff line number Diff line change
Expand Up @@ -18,6 +18,7 @@

#include "../tests.h"
#include "dof_tools_common.h"
#include "dof_tools_common_fake_hp.h"

// check
// DoFTools::
Expand All @@ -26,9 +27,9 @@



template <int dim>
template <typename DoFHandlerType>
void
check_this(const DoFHandler<dim> &dof_handler)
check_this(const DoFHandlerType &dof_handler)
{
// we split up the matrix into
// blocks according to the number
Expand Down
5 changes: 3 additions & 2 deletions tests/dofs/dof_tools_01d.cc
Original file line number Diff line number Diff line change
Expand Up @@ -18,6 +18,7 @@

#include "../tests.h"
#include "dof_tools_common.h"
#include "dof_tools_common_fake_hp.h"

// check
// DoFTools::
Expand All @@ -26,9 +27,9 @@



template <int dim>
template <typename DoFHandlerType>
void
check_this(const DoFHandler<dim> &dof_handler)
check_this(const DoFHandlerType &dof_handler)
{
// we split up the matrix into
// blocks according to the number
Expand Down

0 comments on commit 30eb8db

Please sign in to comment.