-New: DoFTools::extract_dofs_with_support_on() returns a set of degrees of
-freedom that have support on a subset of locally owned cells.
+New: DoFTools::extract_dofs_with_support_contained_within() returns a set of
+degrees of freedom whose support is entirely contained within the cells for
+which the predicate returns true.
<br>
(Denis Davydov, 2017/03/24)
const std::set<types::boundary_id> &boundary_ids = std::set<types::boundary_id>());
/**
- * Extract all degrees of freedom (DoFs) that have support only within cells,
- * for which @p predicate is <code>true</code>. The result is returned as
- * an IndexSet.
+ * Extract all indices of shape functions such that their support is entirely
+ * contained within the cells for which the @p predicate is <code>true</code>.
+ * The result is returned as an IndexSet.
*
* Consider the following FE space where predicate returns <code>true</code>
- * for all cells on the left side:
+ * for all cells on the left half of the domain:
*
- * @image html extract_dofs_with_support_on.png
+ * @image html extract_dofs_with_support_contained_within.png
*
- * This functions will return the union of all DoFs on those cells minus
+ * This functions will return the union of all DoF indices on those cells minus
* DoF 11, 13, 2 and 0; the result will be <code>[9,10], 12, [14,38]</code>.
- * On the image above the returned DoFs are separated from the rest by the
+ * In the image above the returned DoFs are separated from the rest by the
* red line
*
* Essentially, the question this functions answers is the following:
*
* In case of parallel::distributed::Triangulation @p predicate will be called
* only for locally owned and ghost cells. The resulting index set may contain
- * DoFs that are associated to the locally owned or ghost cells, but are not
+ * DoFs that are associated with the locally owned or ghost cells, but are not
* owned by the current MPI core.
*/
template <typename DoFHandlerType>
IndexSet
- extract_dofs_with_support_on(const DoFHandlerType &dof_handler,
- const std::function< bool(const typename DoFHandlerType::active_cell_iterator &)> &predicate,
- const ConstraintMatrix &cm = ConstraintMatrix());
+ extract_dofs_with_support_contained_within(const DoFHandlerType &dof_handler,
+ const std::function< bool(const typename DoFHandlerType::active_cell_iterator &)> &predicate,
+ const ConstraintMatrix &cm = ConstraintMatrix());
/**
* Extract a vector that represents the constant modes of the DoFHandler for
template <typename DoFHandlerType>
IndexSet
- extract_dofs_with_support_on(const DoFHandlerType &dof_handler,
- const std::function< bool(const typename DoFHandlerType::active_cell_iterator &)> &predicate,
- const ConstraintMatrix &cm)
+ extract_dofs_with_support_contained_within(const DoFHandlerType &dof_handler,
+ const std::function< bool(const typename DoFHandlerType::active_cell_iterator &)> &predicate,
+ const ConstraintMatrix &cm)
{
const std::function< bool(const typename DoFHandlerType::active_cell_iterator &)> predicate_local
= [=] (const typename DoFHandlerType::active_cell_iterator & cell) -> bool
}
// Get halo layer and accumulate its DoFs
- std::set<types::global_dof_index> halo_dofs;
+ std::set<types::global_dof_index> dofs_with_support_on_halo_cells;
const std::vector<typename DoFHandlerType::active_cell_iterator> halo_cells =
GridTools::compute_active_cell_halo_layer(dof_handler, predicate_local);
const unsigned int dofs_per_cell = (*it)->get_fe().dofs_per_cell;
local_dof_indices.resize (dofs_per_cell);
(*it)->get_dof_indices (local_dof_indices);
- halo_dofs.insert(local_dof_indices.begin(),
- local_dof_indices.end());
+ dofs_with_support_on_halo_cells.insert(local_dof_indices.begin(),
+ local_dof_indices.end());
}
// A complication coming from the fact that DoFs living on locally
- // owned cells which satisfy predicate may participate in constrains for
+ // owned cells which satisfy predicate may participate in constraints for
// DoFs outside of this region.
if (cm.n_constraints() > 0)
{
- const std::vector<std::pair<types::global_dof_index,double> > *line_ptr;
-
- // Remove DoFs that are part of constrains for DoFs outside
+ // Remove DoFs that are part of constraints for DoFs outside
// of predicate. Since we will subtract halo_dofs from predicate_dofs,
// achieve this by extending halo_dofs with DoFs to which
// halo_dofs are constrained.
std::set<types::global_dof_index> extra_halo;
- for (std::set<types::global_dof_index>::iterator it = halo_dofs.begin();
- it != halo_dofs.end(); ++it)
- {
- line_ptr = cm.get_constraint_entries(*it);
- // if halo DoF is constrained, add all DoFs to which it's constrained.
- if (line_ptr!=NULL)
- {
- const unsigned int line_size = line_ptr->size();
- for (unsigned int j=0; j<line_size; ++j)
- extra_halo.insert((*line_ptr)[j].first);
- }
- }
+ for (std::set<types::global_dof_index>::iterator it = dofs_with_support_on_halo_cells.begin();
+ it != dofs_with_support_on_halo_cells.end(); ++it)
+ // if halo DoF is constrained, add all DoFs to which it's constrained
+ // because after resolving constraints, the support of the DoFs that
+ // constrain the current DoF will extend to the halo cells.
+ if (const std::vector<std::pair<types::global_dof_index,double> > *line_ptr = cm.get_constraint_entries(*it))
+ {
+ const unsigned int line_size = line_ptr->size();
+ for (unsigned int j=0; j<line_size; ++j)
+ extra_halo.insert((*line_ptr)[j].first);
+ }
- halo_dofs.insert(extra_halo.begin(),
- extra_halo.end());
+ dofs_with_support_on_halo_cells.insert(extra_halo.begin(),
+ extra_halo.end());
}
// Rework our std::set's into IndexSet and subtract halo DoFs from the
// predicate_dofs:
IndexSet support_set(dof_handler.n_dofs());
- support_set.add_indices(predicate_dofs.begin(), predicate_dofs.end());
+ support_set.add_indices(predicate_dofs.begin(),
+ predicate_dofs.end());
support_set.compress();
IndexSet halo_set(dof_handler.n_dofs());
- halo_set.add_indices(halo_dofs.begin(), halo_dofs.end());
+ halo_set.add_indices(dofs_with_support_on_halo_cells.begin(),
+ dofs_with_support_on_halo_cells.end());
halo_set.compress();
support_set.subtract_set(halo_set);
// we intentionally do not want to limit the output to locally owned DoFs.
- // In other words, the output may contain DoFs that are on the locally
- // owned cells, but are not owned by this MPI core.
-
return support_set;
}
template
IndexSet
- DoFTools::extract_dofs_with_support_on<DoFHandler<deal_II_dimension> >
+ DoFTools::extract_dofs_with_support_contained_within<DoFHandler<deal_II_dimension> >
(const DoFHandler<deal_II_dimension> &,
const std::function< bool(const typename DoFHandler<deal_II_dimension>::active_cell_iterator &)> &,
const ConstraintMatrix&);
ConstraintMatrix cm;
DoFTools::make_hanging_node_constraints(dh,cm);
- IndexSet support = DoFTools::extract_dofs_with_support_on(dh, pred_d<dim>, cm);
+ IndexSet support = DoFTools::extract_dofs_with_support_contained_within(dh, pred_d<dim>, cm);
support.print(deallog);
// print grid and DoFs for visual inspection
// ---------------------------------------------------------------------
//
-// similar to dof_tools_23, but for the p::d::Tria. As an MPI test, make sure
-// that the total number of shape functions that are non-zero within the domain
-// is the same.
+// test DoFTools::extract_dofs_with_support_contained_within() for the p::d::Tria.
+// As an MPI test, make sure that the total number of shape functions that are
+// non-zero within the domain is the same.
#include "../tests.h"
#include <deal.II/distributed/tria.h>
DoFTools::make_hanging_node_constraints(dh,cm);
cm.close ();
- const IndexSet support = DoFTools::extract_dofs_with_support_on(dh, pred_d<dim>, cm);
+ const IndexSet support = DoFTools::extract_dofs_with_support_contained_within(dh, pred_d<dim>, cm);
const IndexSet support_local = support & dh.locally_owned_dofs();
deallog << support.n_elements() << std::endl;
// ---------------------------------------------------------------------
//
-// same as dof_tools_32 test DoFTools::extract_dofs_with_support_on() but
+// test DoFTools::extract_dofs_with_support_contained_within() but
// for a slightly different configuration
#include "../tests.h"
DoFTools::make_hanging_node_constraints(dh,cm);
IndexSet support = left ?
- DoFTools::extract_dofs_with_support_on(dh, pred_left<dim>, cm) :
- DoFTools::extract_dofs_with_support_on(dh, pred_right<dim>, cm);
+ DoFTools::extract_dofs_with_support_contained_within(dh, pred_left<dim>, cm) :
+ DoFTools::extract_dofs_with_support_contained_within(dh, pred_right<dim>, cm);
support.print(deallog);
// print grid and DoFs for visual inspection
// ---------------------------------------------------------------------
//
-// test DoFTools::extract_dofs_with_support_on() for calculation of the RHS
-// function in parallel
+// test DoFTools::extract_dofs_with_support_contained_within() for calculation
+// of the RHS function in parallel
#include "../tests.h"
#include <deal.II/base/logstream.h>
std::vector<types::global_dof_index> local_dof_indices(fe.dofs_per_cell);
// get support on the predicate
- IndexSet support = DoFTools::extract_dofs_with_support_on(dh, pred_d<dim>, cm);
+ IndexSet support = DoFTools::extract_dofs_with_support_contained_within(dh, pred_d<dim>, cm);
IndexSet local_support = support & locally_owned_set;
// rhs vectors: