Hello deal.II community,
I'm having an issue with what something I believe to be very simple. I have
an existing AffineConstraints "constraints", from which I want to remove a
subset of constraints. This subset is the set of fluid velocity dofs on
some boundary:
IndexSet velocity_dofs =
DoFTools::extract_boundary_dofs(dof_handler, velocity_mask, {id});
IndexSet subset = constraints.get_local_lines();
subset.subtract_set(velocity_dofs);
Then I try to re-create the constraints using get_view() and merge() as
indicated in affine_constraints.h :
auto tmp_constraints = constraints.get_view(subset);
constraints.reinit(locally_owned_dofs, locally_relevant_dofs);
constraints.close();
constraints.merge(tmp_constraints);
But I get this error from get_view(), in debug and on a single MPI proc:
--------------------------------------------------------
An error occurred in line <2589> of file
</home/bawina/Code/dealii/include/deal.II/lac/affine_constraints.h> in
function
dealii::types::global_dof_index
dealii::AffineConstraints<number>::calculate_line_index(size_type) const
[with number = double; dealii::types::global_dof_index = unsigned int;
size_type = unsigned int]
The violated condition was:
local_lines.is_element(line_n)
Additional information:
The index set given to this constraints object indicates constraints
for degree of freedom 10371 should not be stored by this object, but a
constraint is being added.
Stacktrace:
-----------
#0 ./monolithic_fsi:
dealii::AffineConstraints<double>::calculate_line_index(unsigned int) const
#1 ./monolithic_fsi:
dealii::AffineConstraints<double>::is_constrained(unsigned int) const
#2 /usr/local/lib/libdeal_II.g.so.9.7.0-pre:
dealii::AffineConstraints<double>::add_constraint(unsigned int,
dealii::ArrayView<std::pair<unsigned int, double> const,
dealii::MemorySpace::Host> const&, double)
#3 /usr/local/lib/libdeal_II.g.so.9.7.0-pre:
dealii::AffineConstraints<double>::get_view(dealii::IndexSet const&) const
[...]
It comes from view.add_constraint(mask.index_within_set(line.index), ...,
...) in get_view(), is it a "shift" issue?
Copying and filtering by hand works:
{
AffineConstraints<double> filtered;
filtered.reinit(locally_owned_dofs, locally_relevant_dofs);
for (const auto &line : constraints.get_lines())
{
if (velocity_dofs.is_element(line.index))
continue;
filtered.add_line(line.index);
filtered.add_entries(line.index, line.entries);
filtered.set_inhomogeneity(line.index, line.inhomogeneity);
for (const auto &entry : line.entries)
AssertThrow(!velocity_dofs.is_element(entry.first),
ExcMessage("Constraint involve a boundary velocity
dof"));
}
filtered.close();
constraints = std::move(filtered);
}
Thank you for your time!
Arthur
--
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see
https://groups.google.com/d/forum/dealii?hl=en
---
You received this message because you are subscribed to the Google Groups
"deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email
to [email protected].
To view this discussion visit
https://groups.google.com/d/msgid/dealii/f1b66d74-5768-42d5-9f47-45a07a1b9cd9n%40googlegroups.com.