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

Iterator filter predicates for boundary and manifold IDs #13154

Merged
merged 4 commits into from
Dec 30, 2021
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
5 changes: 5 additions & 0 deletions doc/news/changes/minor/20211230Jean-PaulPelteret
Original file line number Diff line number Diff line change
@@ -0,0 +1,5 @@
New: The IteratorFilters::BoundaryIdEqualTo and IteratorFilters::ManifoldIdEqualTo
filtered iterator predicates have been added to help select cell faces that
have, respectively, one or more specified boundary ID(s) and manifold ID(s).
<br>
(Jean-Paul Pelteret, 2021/12/30)
124 changes: 124 additions & 0 deletions include/deal.II/grid/filtered_iterator.h
Original file line number Diff line number Diff line change
Expand Up @@ -328,6 +328,82 @@ namespace IteratorFilters
bool
operator()(const Iterator &i) const;
};


/**
* Filter for iterators that evaluates to true if the iterator of the object
* pointed to is equal to a value or set of values given to the constructor,
* assuming that the iterator allows querying for a boundary id.
*
* @ingroup Iterators
*/
class BoundaryIdEqualTo
{
public:
/**
* Constructor. Store the boundary id which iterators shall have to be
* evaluated to true.
*/
BoundaryIdEqualTo(const types::boundary_id boundary_id);

/**
* Constructor. Store a collection of boundary ids which iterators shall
* have to be evaluated to true.
*/
BoundaryIdEqualTo(const std::set<types::boundary_id> &boundary_ids);

/**
* Evaluation operator. Returns true if the boundary id of the object
* pointed to is equal within the stored set of value allowable values.
*/
template <class Iterator>
bool
operator()(const Iterator &i) const;

protected:
/**
* Stored value to compare the material id with.
*/
const std::set<types::boundary_id> boundary_ids;
};


/**
* Filter for iterators that evaluates to true if the iterator of the object
* pointed to is equal to a value or set of values given to the constructor,
* assuming that the iterator allows querying for a manifold id.
*
* @ingroup Iterators
*/
class ManifoldIdEqualTo
{
public:
/**
* Constructor. Store the boundary id which iterators shall have to be
* evaluated to true.
*/
ManifoldIdEqualTo(const types::manifold_id manifold_id);

/**
* Constructor. Store a collection of boundary ids which iterators shall
* have to be evaluated to true.
*/
ManifoldIdEqualTo(const std::set<types::manifold_id> &manifold_ids);

/**
* Evaluation operator. Returns true if the boundary id of the object
* pointed to is equal within the stored set of value allowable values.
*/
template <class Iterator>
bool
operator()(const Iterator &i) const;

protected:
/**
* Stored value to compare the material id with.
*/
const std::set<types::manifold_id> manifold_ids;
};
} // namespace IteratorFilters


Expand Down Expand Up @@ -1440,6 +1516,54 @@ namespace IteratorFilters
{
return (i->at_boundary());
}



// ---------------- IteratorFilters::BoundaryIdEqualTo ---------
inline BoundaryIdEqualTo::BoundaryIdEqualTo(
const types::boundary_id boundary_id)
: boundary_ids{boundary_id}
{}



inline BoundaryIdEqualTo::BoundaryIdEqualTo(
const std::set<types::boundary_id> &boundary_ids)
: boundary_ids(boundary_ids)
{}



template <class Iterator>
inline bool
BoundaryIdEqualTo::operator()(const Iterator &i) const
{
return boundary_ids.find(i->boundary_id()) != boundary_ids.end();
}



// ---------------- IteratorFilters::ManifoldIdEqualTo ---------
inline ManifoldIdEqualTo::ManifoldIdEqualTo(
const types::manifold_id manifold_id)
: manifold_ids{manifold_id}
{}



inline ManifoldIdEqualTo::ManifoldIdEqualTo(
const std::set<types::manifold_id> &manifold_ids)
: manifold_ids(manifold_ids)
{}



template <class Iterator>
inline bool
ManifoldIdEqualTo::operator()(const Iterator &i) const
{
return manifold_ids.find(i->manifold_id()) != manifold_ids.end();
}
} // namespace IteratorFilters


Expand Down
53 changes: 37 additions & 16 deletions tests/grid/filtered_iterator_06.cc
Original file line number Diff line number Diff line change
Expand Up @@ -31,7 +31,8 @@ void
test()
{
Triangulation<dim, spacedim> tria;
GridGenerator::hyper_cube(tria);
GridGenerator::hyper_cube(tria, 0, 1, true);

// refine the boundary cells a few times
for (unsigned int i = 0; i < 5; ++i)
{
Expand All @@ -48,29 +49,49 @@ test()
++i;
}

// Count the faces that are on the boundary and have a manifold_id of 0
// Count the faces that are on the boundary and have a boundary_id of 0 and a
// manifold_id of 0
const types::boundary_id boundary_id = 0;
const types::manifold_id manifold_id = 0;
std::set<typename Triangulation<dim, spacedim>::active_face_iterator>
face_set;
boundary_face_set, manifold_face_set;

for (const auto &cell : tria.active_cell_iterators())
for (const unsigned int face_n : GeometryInfo<dim>::face_indices())
if (cell->face(face_n)->at_boundary() &&
cell->face(face_n)->manifold_id() == manifold_id)
face_set.insert(cell->face(face_n));
if (cell->face(face_n)->at_boundary())
{
if (cell->face(face_n)->boundary_id() == boundary_id)
boundary_face_set.insert(cell->face(face_n));

if (cell->face(face_n)->manifold_id() == manifold_id)
manifold_face_set.insert(cell->face(face_n));
}

std::size_t n_boundary_filtered_cells = 0;
for (const auto &filtered_cell :
filter_iterators(tria.active_face_iterators(),
IteratorFilters::AtBoundary(),
IteratorFilters::BoundaryIdEqualTo(boundary_id)))
{
AssertThrow(boundary_face_set.count(filtered_cell) == 1,
ExcMessage("Wrong cell filtered."));
++n_boundary_filtered_cells;
}
AssertThrow(n_boundary_filtered_cells == boundary_face_set.size(),
ExcMessage("Boundary filtered cells missing."));

std::size_t n_filtered_cells = 0;
for (const auto filtered_cell : filter_iterators(
tria.active_face_iterators(),
IteratorFilters::AtBoundary(),
[=](const typename Triangulation<dim, spacedim>::active_face_iterator
&face) { return face->manifold_id() == manifold_id; }))
std::size_t n_manifold_filtered_cells = 0;
for (const auto &filtered_cell :
filter_iterators(tria.active_face_iterators(),
IteratorFilters::AtBoundary(),
IteratorFilters::ManifoldIdEqualTo(manifold_id)))
{
AssertThrow(face_set.count(filtered_cell) == 1,
AssertThrow(manifold_face_set.count(filtered_cell) == 1,
ExcMessage("Wrong cell filtered."));
++n_filtered_cells;
++n_manifold_filtered_cells;
}
AssertThrow(n_filtered_cells == face_set.size(),
ExcMessage("Filtered cells missing."));
AssertThrow(n_manifold_filtered_cells == manifold_face_set.size(),
ExcMessage("Manifold filtered cells missing."));
}

int
Expand Down
48 changes: 34 additions & 14 deletions tests/grid/filtered_iterator_06_operator.cc
Original file line number Diff line number Diff line change
Expand Up @@ -34,7 +34,8 @@ void
test()
{
Triangulation<dim, spacedim> tria;
GridGenerator::hyper_cube(tria);
GridGenerator::hyper_cube(tria, 0, 1, true);

// refine the boundary cells a few times
for (unsigned int i = 0; i < 5; ++i)
{
Expand All @@ -51,28 +52,47 @@ test()
++i;
}

// Count the faces that are on the boundary and have a manifold_id of 0
// Count the faces that are on the boundary and have a boundary_id of 0 and a
// manifold_id of 0
const types::boundary_id boundary_id = 0;
const types::manifold_id manifold_id = 0;
std::set<typename Triangulation<dim, spacedim>::active_face_iterator>
face_set;
boundary_face_set, manifold_face_set;

for (const auto &cell : tria.active_cell_iterators())
for (const unsigned int face_n : GeometryInfo<dim>::face_indices())
if (cell->face(face_n)->at_boundary() &&
cell->face(face_n)->manifold_id() == manifold_id)
face_set.insert(cell->face(face_n));
if (cell->face(face_n)->at_boundary())
{
if (cell->face(face_n)->boundary_id() == boundary_id)
boundary_face_set.insert(cell->face(face_n));

if (cell->face(face_n)->manifold_id() == manifold_id)
manifold_face_set.insert(cell->face(face_n));
}

std::size_t n_boundary_filtered_cells = 0;
for (const auto &filtered_cell :
tria.active_face_iterators() | IteratorFilters::AtBoundary() |
IteratorFilters::BoundaryIdEqualTo(boundary_id))
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

:-) That's just elegant!

Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Yeah, I'm really happy that you implemented operator|. It's so nice to chain these filters!

Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I really like that syntax. It's been used on the unix command line for 50 years, of course, and I like that we'll be able to use it in C++ as well (or at least will be with C++20).

{
AssertThrow(boundary_face_set.count(filtered_cell) == 1,
ExcMessage("Wrong cell filtered."));
++n_boundary_filtered_cells;
}
AssertThrow(n_boundary_filtered_cells == boundary_face_set.size(),
ExcMessage("Boundary filtered cells missing."));

std::size_t n_filtered_cells = 0;
for (const auto filtered_cell :
std::size_t n_manifold_filtered_cells = 0;
for (const auto &filtered_cell :
tria.active_face_iterators() | IteratorFilters::AtBoundary() |
[=](const typename Triangulation<dim, spacedim>::active_face_iterator
&face) { return face->manifold_id() == manifold_id; })
IteratorFilters::ManifoldIdEqualTo(manifold_id))
{
AssertThrow(face_set.count(filtered_cell) == 1,
AssertThrow(manifold_face_set.count(filtered_cell) == 1,
ExcMessage("Wrong cell filtered."));
++n_filtered_cells;
++n_manifold_filtered_cells;
}
AssertThrow(n_filtered_cells == face_set.size(),
ExcMessage("Filtered cells missing."));
AssertThrow(n_manifold_filtered_cells == manifold_face_set.size(),
ExcMessage("Manifold filtered cells missing."));
}

int
Expand Down