static constexpr std::array<int, 2> unit_normal_orientation{{-1, 1}};
+ static constexpr std::array<Tensor<1, 1>, 2> unit_normal_vector{
+ {Tensor<1, 1>{{-1}}, Tensor<1, 1>{{1}}}};
+
static constexpr std::array<unsigned int, 2> opposite_face{{1, 0}};
static constexpr std::array<unsigned int, 2> dx_to_deal{{0, 1}};
static constexpr std::array<int, 4> unit_normal_orientation{
{-1, 1, -1, 1}};
+ static constexpr std::array<Tensor<1, 2>, 4> unit_normal_vector{
+ {Tensor<1, 2>{{-1., 0.}},
+ Tensor<1, 2>{{1., 0.}},
+ Tensor<1, 2>{{0., -1.}},
+ Tensor<1, 2>{{0., 1.}}}};
+
static constexpr std::array<unsigned int, 4> opposite_face{{1, 0, 3, 2}};
static constexpr std::array<unsigned int, 4> dx_to_deal{{0, 2, 1, 3}};
static constexpr std::array<int, 6> unit_normal_orientation{
{-1, 1, -1, 1, -1, 1}};
+ static constexpr std::array<Tensor<1, 3>, 6> unit_normal_vector{
+ {Tensor<1, 3>{{-1, 0, 0}},
+ Tensor<1, 3>{{1, 0, 0}},
+ Tensor<1, 3>{{0, -1, 0}},
+ Tensor<1, 3>{{0, 1, 0}},
+ Tensor<1, 3>{{0, 0, -1}},
+ Tensor<1, 3>{{0, 0, 1}}}};
+
static constexpr std::array<unsigned int, 6> opposite_face{
{1, 0, 3, 2, 5, 4}};
static constexpr std::array<int, 8> unit_normal_orientation{
{-1, 1, -1, 1, -1, 1, -1, 1}};
+ static constexpr std::array<Tensor<1, 4>, 8> unit_normal_vector{
+ {Tensor<1, 4>{{-1, 0, 0, 0}},
+ Tensor<1, 4>{{1, 0, 0, 0}},
+ Tensor<1, 4>{{0, -1, 0, 0}},
+ Tensor<1, 4>{{0, 1, 0, 0}},
+ Tensor<1, 4>{{0, 0, -1, 0}},
+ Tensor<1, 4>{{0, 0, 1, 0}},
+ Tensor<1, 4>{{0, 0, 0, -1}},
+ Tensor<1, 4>{{0, 0, 0, 1}}}};
+
static constexpr std::array<unsigned int, 8> opposite_face{
{1, 0, 3, 2, 5, 4, 7, 6}};
static constexpr std::array<int, faces_per_cell> unit_normal_orientation =
internal::GeometryInfoHelper::Initializers<dim>::unit_normal_orientation;
+ /**
+ * Unit normal vector (Point<dim>) of a face of the reference cell.
+ *
+ * Note that this is only the <em>standard orientation</em> of faces. At
+ * least in 3d, actual faces of cells in a triangulation can also have the
+ * opposite orientation, depending on a flag that one can query from the
+ * cell it belongs to. For more information, see the
+ * @ref GlossFaceOrientation "glossary"
+ * entry on face orientation.
+ */
+ static constexpr std::array<Tensor<1, dim>, faces_per_cell>
+ unit_normal_vector =
+ internal::GeometryInfoHelper::Initializers<dim>::unit_normal_vector;
+
/**
* List of numbers which denotes which face is opposite to a given face. Its
* entries are the first <tt>2*dim</tt> entries of <tt>{ 1, 0, 3, 2, 5, 4,
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2019 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 GeometryInfo<dim>::unit_normal_vector
+
+#include <deal.II/base/geometry_info.h>
+
+#include "../tests.h"
+
+const double eps = 1e-10;
+DeclException3(ExcWrongValue,
+ double,
+ double,
+ double,
+ << arg1 << " != " << arg2 << " with delta = " << arg3);
+
+template <int dim>
+void
+test()
+{
+ for (unsigned int i = 0; i < GeometryInfo<dim>::faces_per_cell; ++i)
+ {
+ const Tensor<1, dim> unit_normal_vector =
+ GeometryInfo<dim>::unit_normal_vector[i];
+ deallog << "Direction " << i << " = " << unit_normal_vector << std::endl;
+ // check consistency with other arrays within GeometryInfo
+ for (unsigned int j = 0; j < dim; j++)
+ {
+ const double obtained = unit_normal_vector[j];
+ const double expected =
+ (j == GeometryInfo<dim>::unit_normal_direction[i]) ?
+ static_cast<double>(
+ GeometryInfo<dim>::unit_normal_orientation[i]) :
+ 0;
+ AssertThrow(obtained == expected,
+ ExcWrongValue(obtained, expected, obtained - expected));
+ }
+ }
+ deallog << "OK" << std::endl;
+}
+
+
+int
+main()
+{
+ initlog();
+
+ deallog.push("1d");
+ test<1>();
+ deallog.pop();
+ deallog.push("2d");
+ test<2>();
+ deallog.pop();
+ deallog.push("3d");
+ test<3>();
+ deallog.pop();
+ return 0;
+}