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

Fix face orientation #15678

Merged
merged 1 commit into from
Jul 8, 2023
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
12 changes: 6 additions & 6 deletions include/deal.II/grid/reference_cell.h
Original file line number Diff line number Diff line change
Expand Up @@ -875,9 +875,9 @@ class ReferenceCell
{{{0, 2, 1}},
{{0, 1, 2}},
{{2, 1, 0}},
{{1, 2, 0}},
{{2, 0, 1}},
{{1, 0, 2}},
{{2, 0, 1}}}};
{{1, 2, 0}}}};

/**
* Table containing all vertex permutations for a quadrilateral.
Expand Down Expand Up @@ -2016,9 +2016,9 @@ ReferenceCell::standard_to_real_face_line(
static constexpr ndarray<unsigned int, 6, 3> triangle_table = {{{{2, 1, 0}},
{{0, 1, 2}},
{{1, 0, 2}},
{{1, 2, 0}},
{{2, 0, 1}},
{{0, 2, 1}},
{{2, 0, 1}}}};
{{1, 2, 0}}}};


switch (this->kind)
Expand Down Expand Up @@ -2750,11 +2750,11 @@ ReferenceCell::get_combined_orientation(

// face_orientation=true, face_rotation=true, face_flip=false
if (v0_equals({vertices_1[1], vertices_1[2], vertices_1[0]}))
return 3;
return 5;

// face_orientation=true, face_rotation=false, face_flip=true
if (v0_equals({vertices_1[2], vertices_1[0], vertices_1[1]}))
return 5;
return 3;

// face_orientation=false, face_rotation=false, face_flip=false
if (v0_equals({vertices_1[0], vertices_1[2], vertices_1[1]}))
Expand Down
20 changes: 10 additions & 10 deletions tests/simplex/orientation_02.cc
Original file line number Diff line number Diff line change
Expand Up @@ -47,11 +47,8 @@ test(const unsigned int orientation)
const auto &face = dummy.begin()->face(face_no);
const auto permuted =
ReferenceCell(ReferenceCells::Triangle)
.permute_according_orientation(
std::array<unsigned int, 3>{{face->vertex_index(0),
face->vertex_index(1),
face->vertex_index(2)}},
orientation);
.permute_according_orientation(std::array<unsigned int, 3>{{0, 1, 2}},
orientation);

auto direction =
cross_product_3d(vertices[permuted[1]] - vertices[permuted[0]],
Expand All @@ -61,7 +58,12 @@ test(const unsigned int orientation)
vertices.push_back(face->center() + direction);

CellData<3> cell;
cell.vertices = {permuted[0], permuted[1], permuted[2], 4u};
cell.vertices.resize(4);

cell.vertices[permuted[0]] = face->vertex_index(0);
cell.vertices[permuted[1]] = face->vertex_index(1);
cell.vertices[permuted[2]] = face->vertex_index(2);
cell.vertices[3] = 4;

cell.material_id = 1;
cells.push_back(cell);
Expand All @@ -73,10 +75,8 @@ test(const unsigned int orientation)
cell++;

// check orientation
deallog << "face orientation: " << orientation << ' ' << std::endl;
AssertDimension(orientation,
(cell->face_orientation(0) * 1 + cell->face_rotation(0) * 2 +
cell->face_flip(0) * 4));
deallog << "face orientation: " << orientation << ' '
<< int(cell->combined_face_orientation(0)) << ' ' << std::endl;

// check vertices
{
Expand Down
20 changes: 10 additions & 10 deletions tests/simplex/orientation_02.output
Original file line number Diff line number Diff line change
@@ -1,37 +1,37 @@

DEAL::face orientation: 0
DEAL::face orientation: 0 0
DEAL::0 2 1 4 vs. 0 2 1 4
DEAL::0-2 vs. 0-2
DEAL::1-2 vs. 1-2
DEAL::0-1 vs. 0-1
DEAL::
DEAL::face orientation: 1
DEAL::face orientation: 1 1
DEAL::0 1 2 4 vs. 0 1 2 4
DEAL::0-1 vs. 0-1
DEAL::1-2 vs. 1-2
DEAL::0-2 vs. 0-2
DEAL::
DEAL::face orientation: 2
DEAL::face orientation: 2 2
DEAL::2 1 0 4 vs. 2 1 0 4
DEAL::1-2 vs. 1-2
DEAL::0-1 vs. 0-1
DEAL::0-2 vs. 0-2
DEAL::
DEAL::face orientation: 3
DEAL::1 2 0 4 vs. 1 2 0 4
DEAL::1-2 vs. 1-2
DEAL::face orientation: 3 3
DEAL::2 0 1 4 vs. 2 0 1 4
DEAL::0-2 vs. 0-2
DEAL::0-1 vs. 0-1
DEAL::1-2 vs. 1-2
DEAL::
DEAL::face orientation: 4
DEAL::face orientation: 4 4
DEAL::1 0 2 4 vs. 1 0 2 4
DEAL::0-1 vs. 0-1
DEAL::0-2 vs. 0-2
DEAL::1-2 vs. 1-2
DEAL::
DEAL::face orientation: 5
DEAL::2 0 1 4 vs. 2 0 1 4
DEAL::face orientation: 5 5
DEAL::1 2 0 4 vs. 1 2 0 4
DEAL::1-2 vs. 1-2
DEAL::0-2 vs. 0-2
DEAL::0-1 vs. 0-1
DEAL::1-2 vs. 1-2
DEAL::
175 changes: 175 additions & 0 deletions tests/simplex/orientation_03.cc
Original file line number Diff line number Diff line change
@@ -0,0 +1,175 @@
// ---------------------------------------------------------------------
//
// Copyright (C) 2023 by the deal.II authors
//
// This file is part of the deal.II library.
//
// The deal.II library is free software; you can use it, redistribute
// it, and/or modify it under the terms of the GNU Lesser General
// Public License as published by the Free Software Foundation; either
// version 2.1 of the License, or (at your option) any later version.
// The full text of the license can be found in the file LICENSE.md at
// the top level directory of deal.II.
//
// ---------------------------------------------------------------------

// Test a mesh with two tetrahedra for all possible orientations. Similar to
// orientation_02 but also checks that quadrature points on faces (computed via
// FEFaceValues) are correct.
Comment on lines +16 to +18
Copy link
Member

Choose a reason for hiding this comment

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

Nice!


#include <deal.II/fe/fe_nothing.h>
#include <deal.II/fe/fe_simplex_p.h>
#include <deal.II/fe/fe_values.h>
#include <deal.II/fe/mapping.h>

#include <deal.II/grid/grid_generator.h>
#include <deal.II/grid/reference_cell.h>
#include <deal.II/grid/tria.h>

#include "../tests.h"

void
test_face(const std::vector<Point<3>> & vertices_,
const std::vector<CellData<3>> &cell_data_,
const unsigned int face_n)
{
const ReferenceCell ref_cell = ReferenceCells::Tetrahedron;
auto vertices = vertices_;
auto cell_data = cell_data_;

Point<3> extra_vertex;
for (unsigned int i = 0; i < 3; ++i)
extra_vertex += ref_cell.template vertex<3>(ref_cell.face_to_cell_vertices(
face_n, i, ReferenceCell::default_combined_face_orientation()));

extra_vertex /= 3.0;
extra_vertex += ref_cell.template unit_normal_vectors<3>(face_n);

vertices.push_back(extra_vertex);

cell_data.emplace_back();
cell_data.back().vertices.resize(0);
for (unsigned int i = 0; i < 3; ++i)
cell_data.back().vertices.push_back(ref_cell.face_to_cell_vertices(
face_n, i, ref_cell.default_combined_face_orientation()));
cell_data.back().vertices.push_back(ref_cell.n_vertices());
std::sort(cell_data.back().vertices.begin(), cell_data.back().vertices.end());

unsigned int permutation_n = 0;
do
{
Triangulation<3> tria;
tria.create_triangulation(vertices, cell_data, SubCellData());
deallog << " p2 = " << permutation_n;

FE_Nothing<3> fe(ref_cell);
const Mapping<3> &mapping =
ref_cell.template get_default_linear_mapping<3>();
Quadrature<2> face_quad(FE_SimplexP<3>(1).get_unit_face_support_points());

FEFaceValues<3> cell_face_values(mapping,
fe,
face_quad,
update_quadrature_points);
FEFaceValues<3> neighbor_cell_face_values(mapping,
fe,
face_quad,
update_quadrature_points);

for (const auto &cell : tria.active_cell_iterators())
for (unsigned int face_no : cell->face_indices())
{
auto neighbor_cell = cell->neighbor(face_no);
if (neighbor_cell == tria.end())
continue;

auto face = cell->face(face_no);
cell_face_values.reinit(cell, face);
unsigned int neighbor_face_no = 0;
for (; neighbor_face_no < ref_cell.n_faces(); ++neighbor_face_no)
if (neighbor_cell->face(neighbor_face_no) == face)
break;
neighbor_cell_face_values.reinit(neighbor_cell, neighbor_face_no);

deallog << " : " << int(cell->combined_face_orientation(face_no))
<< ", "
<< int(neighbor_cell->combined_face_orientation(
neighbor_face_no));

double max_distance = 0.0;
for (unsigned int q = 0; q < face_quad.size(); ++q)
max_distance = std::max(
max_distance,
cell_face_values.get_quadrature_points()[q].distance(
neighbor_cell_face_values.get_quadrature_points()[q]));

if (max_distance > 1e-12)
{
deallog << "!!!!! Found a wrong point permutation !!!!!"
<< std::endl;

deallog << "cell points =" << std::endl;
for (unsigned int q = 0; q < face_quad.size(); ++q)
{
deallog << " "
<< cell_face_values.get_quadrature_points()[q]
<< std::endl;
}
deallog << std::endl;

deallog << "neighbor_cell points =" << std::endl;
for (unsigned int q = 0; q < face_quad.size(); ++q)
{
deallog
<< " "
<< neighbor_cell_face_values.get_quadrature_points()[q]
<< std::endl;
}
deallog << std::endl;

AssertThrow(false, ExcMessage("Should match"));
}
}
deallog << std::endl;
++permutation_n;
}
while (std::next_permutation(cell_data.back().vertices.begin(),
cell_data.back().vertices.end()));
}

void
test()
{
std::vector<Point<3>> vertices;
vertices.emplace_back(0.0, 0.0, 0.0);
vertices.emplace_back(1.0, 0.0, 0.0);
vertices.emplace_back(0.0, 1.0, 0.0);
vertices.emplace_back(0.0, 0.0, 1.0);

std::vector<CellData<3>> cells;
cells.emplace_back();
cells.back().vertices = {0, 1, 2, 3};

for (unsigned int face_n = 0; face_n < 4; ++face_n)
{
std::sort(cells.back().vertices.begin(), cells.back().vertices.end());
unsigned int permutation_n = 0;
deallog << "face_n = " << face_n << std::endl;
do
{
deallog << "p1 = " << permutation_n << std::endl;
test_face(vertices, cells, face_n);
++permutation_n;
}
while (std::next_permutation(cells.back().vertices.begin(),
cells.back().vertices.end()));
}
}

int
main()
{
initlog();

test();
}