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

Add two asserts in NonMatching::MeshClassifier. #13709

Merged
merged 2 commits into from
May 10, 2022
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
19 changes: 19 additions & 0 deletions source/non_matching/mesh_classifier.cc
Original file line number Diff line number Diff line change
Expand Up @@ -46,6 +46,16 @@ namespace NonMatching
{
namespace MeshClassifierImplementation
{
DeclExceptionMsg(
ReclassifyNotCalled,
"The Triangulation has not been classified. You need to call the "
"reclassify()-function before using this function.");

DeclExceptionMsg(
ExcTriangulationMismatch,
"The incoming cell does not belong to the triangulation passed to "
"the constructor.");

/**
* Return LocationToLevelSet::inside/outside if all values in incoming
* vector are negative/positive, otherwise return
Expand Down Expand Up @@ -419,6 +429,11 @@ namespace NonMatching
MeshClassifier<dim>::location_to_level_set(
const typename Triangulation<dim>::cell_iterator &cell) const
{
Assert(cell_locations.size() == triangulation->n_active_cells(),
internal::MeshClassifierImplementation::ReclassifyNotCalled());
Assert(&cell->get_triangulation() == triangulation,
internal::MeshClassifierImplementation::ExcTriangulationMismatch());

return cell_locations.at(cell->active_cell_index());
}

Expand All @@ -431,6 +446,10 @@ namespace NonMatching
const unsigned int face_index) const
{
AssertIndexRange(face_index, GeometryInfo<dim>::faces_per_cell);
Assert(face_locations.size() == triangulation->n_raw_faces(),
internal::MeshClassifierImplementation::ReclassifyNotCalled());
Assert(&cell->get_triangulation() == triangulation,
internal::MeshClassifierImplementation::ExcTriangulationMismatch());

return face_locations.at(cell->face(face_index)->index());
}
Expand Down
12 changes: 6 additions & 6 deletions tests/non_matching/mesh_classifier.cc
Original file line number Diff line number Diff line change
Expand Up @@ -73,15 +73,15 @@ print_cell_and_face_locations(
const NonMatching::MeshClassifier<dim> & classifier,
const typename Triangulation<dim>::active_cell_iterator &cell)
{
const NonMatching::LocationToLevelSet cell_position =
const NonMatching::LocationToLevelSet cell_location =
classifier.location_to_level_set(cell);
deallog << "cell " << location_to_string(cell_position) << std::endl;
deallog << "cell " << location_to_string(cell_location) << std::endl;

for (unsigned int i = 0; i < GeometryInfo<dim>::faces_per_cell; ++i)
for (unsigned int f : cell->face_indices())
{
const NonMatching::LocationToLevelSet cell_position =
classifier.location_to_level_set(cell, i);
deallog << "face " << i << ' ' << location_to_string(cell_position)
const NonMatching::LocationToLevelSet face_location =
classifier.location_to_level_set(cell, f);
deallog << "face " << f << ' ' << location_to_string(face_location)
<< std::endl;
}
}
Expand Down