const bool with_constraints = true) const;
/**
- * This internal method takes the local indices on a cell and fills them
- * into this class. It resolves the constraints and distributes the
- * results. Ghost indices, i.e., indices that are located on another
- * processor, get a temporary number by this function, and will later be
- * assigned the correct index after all the ghost indices have been
- * collected by the call to @p assign_ghosts.
+ * This internal method takes the local indices on a cell (two versions:
+ * hanging-node constraints resolved if possible and plain, i.e., not
+ * resolved) and fills them into this class. It resolves the constraints
+ * and distributes the results. Ghost indices, i.e., indices that are
+ * located on another processor, get a temporary number by this function,
+ * and will later be assigned the correct index after all the ghost
+ * indices have been collected by the call to @p assign_ghosts.
*/
template <typename number>
void
read_dof_indices(
+ const std::vector<types::global_dof_index> &local_indices_resolved,
const std::vector<types::global_dof_index> &local_indices,
- const std::vector<unsigned int> & lexicographic_inv,
const dealii::AffineConstraints<number> & constraints,
const unsigned int cell_number,
ConstraintValues<double> & constraint_values,
*/
std::vector<unsigned int> dof_indices;
+ /**
+ * Masks indicating for each cell and component if the optimized
+ * hanging-node constraint is applicable and if yes which type.
+ */
+ std::vector<unsigned int> hanging_node_constraint_masks;
+
/**
* This variable describes the position of constraints in terms of the
* local numbering of degrees of freedom on a cell. The first number
template <typename number>
void
DoFInfo::read_dof_indices(
+ const std::vector<types::global_dof_index> &local_indices_resolved,
const std::vector<types::global_dof_index> &local_indices,
- const std::vector<unsigned int> & lexicographic_inv,
const dealii::AffineConstraints<number> & constraints,
const unsigned int cell_number,
ConstraintValues<double> & constraint_values,
i < component_dof_indices_offset[fe_index][comp + 1];
i++)
{
- types::global_dof_index current_dof =
- local_indices[lexicographic_inv[i]];
- const auto *entries_ptr =
+ types::global_dof_index current_dof = local_indices_resolved[i];
+ const auto * entries_ptr =
constraints.get_constraint_entries(current_dof);
// dof is constrained
{
for (unsigned int i = 0; i < dofs_this_cell; ++i)
{
- types::global_dof_index current_dof =
- local_indices[lexicographic_inv[i]];
+ types::global_dof_index current_dof = local_indices[i];
if (n_mpi_procs > 1 &&
(current_dof < first_owned || current_dof >= last_owned))
{
AssertDimension(n_dof_handlers, locally_owned_dofs.size());
AssertDimension(n_dof_handlers, constraint.size());
- std::vector<types::global_dof_index> local_dof_indices;
+ std::vector<types::global_dof_index> local_dof_indices_resolved;
+ std::vector<types::global_dof_index> local_dof_indices;
std::vector<std::vector<std::vector<unsigned int>>> lexicographic(
n_dof_handlers);
0;
if (dofh->get_fe_collection().size() > 1)
dof_info[no].cell_active_fe_index[counter] = fe_index;
- local_dof_indices.resize(dof_info[no].dofs_per_cell[fe_index]);
- cell_it->get_dof_indices(local_dof_indices);
+ local_dof_indices_resolved.resize(
+ dof_info[no].dofs_per_cell[fe_index]);
+ cell_it->get_dof_indices(local_dof_indices_resolved);
+
+ local_dof_indices.resize(local_dof_indices_resolved.size());
+ for (unsigned int i = 0; i < local_dof_indices_resolved.size();
+ ++i)
+ local_dof_indices[i] =
+ local_dof_indices_resolved[lexicographic[no][fe_index][i]];
+
dof_info[no].read_dof_indices(local_dof_indices,
- lexicographic[no][fe_index],
+ local_dof_indices,
*constraint[no],
counter,
constraint_values,
cell_level_index[counter].first,
cell_level_index[counter].second,
dofh);
- local_dof_indices.resize(dof_info[no].dofs_per_cell[0]);
- cell_it->get_mg_dof_indices(local_dof_indices);
+ local_dof_indices_resolved.resize(
+ dof_info[no].dofs_per_cell[0]);
+ cell_it->get_mg_dof_indices(local_dof_indices_resolved);
+
+ local_dof_indices.resize(local_dof_indices_resolved.size());
+ for (unsigned int i = 0; i < local_dof_indices_resolved.size();
+ ++i)
+ local_dof_indices[i] =
+ local_dof_indices_resolved[lexicographic[no][0][i]];
+
dof_info[no].read_dof_indices(local_dof_indices,
- lexicographic[no][0],
+ local_dof_indices,
*constraint[no],
counter,
constraint_values,
template void
DoFInfo::read_dof_indices<double>(
const std::vector<types::global_dof_index> &,
- const std::vector<unsigned int> &,
+ const std::vector<types::global_dof_index> &,
const dealii::AffineConstraints<double> &,
const unsigned int,
ConstraintValues<double> &,
template void
DoFInfo::read_dof_indices<float>(
const std::vector<types::global_dof_index> &,
- const std::vector<unsigned int> &,
+ const std::vector<types::global_dof_index> &,
const dealii::AffineConstraints<float> &,
const unsigned int,
ConstraintValues<double> &,