const unsigned int (&vertex_indices)[GeometryInfo<dim>::vertices_per_cell]);
/**
- * A version of the last function that can accept input for nonzero
+ * A variant of cell_measure() accepting an ArrayView instead of a
+ * fixed-sized array for @p vertex_indices.
+ *
+ * The parameter @p vertex_indices is expected to have
+ * GeometryInfo<dim>::vertices_per_cell entries. A std::vector is implicitly
+ * convertible to an ArrayView, so it can be passed directly. See the
+ * ArrayView class for more information.
+ */
+ template <int dim>
+ double
+ cell_measure(const std::vector<Point<dim>> & all_vertices,
+ const ArrayView<const unsigned int> &vertex_indices);
+
+ /**
+ * A version of the function above that can accept input for nonzero
* codimension cases. This function only exists to aid generic programming
* and calling it will just raise an exception.
*/
namespace GridTools
{
- // declare specializations
- template <>
- double
- cell_measure<1>(const std::vector<Point<1>> &,
- const unsigned int (&)[GeometryInfo<1>::vertices_per_cell]);
-
- template <>
- double
- cell_measure<2>(const std::vector<Point<2>> &,
- const unsigned int (&)[GeometryInfo<2>::vertices_per_cell]);
-
- template <>
+ template <int dim>
double
- cell_measure<3>(const std::vector<Point<3>> &,
- const unsigned int (&)[GeometryInfo<3>::vertices_per_cell]);
-
-
+ cell_measure(
+ const std::vector<Point<dim>> &all_vertices,
+ const unsigned int (&indices)[GeometryInfo<dim>::vertices_per_cell])
+ {
+ // We forward call to the ArrayView version:
+ const ArrayView<const unsigned int> view(
+ &indices[0], GeometryInfo<dim>::vertices_per_cell);
+ return cell_measure(all_vertices, view);
+ }
template <int dim, typename T>
double
{
template <>
double
- cell_measure<1>(
- const std::vector<Point<1>> &all_vertices,
- const unsigned int (&vertex_indices)[GeometryInfo<1>::vertices_per_cell])
+ cell_measure<1>(const std::vector<Point<1>> & all_vertices,
+ const ArrayView<const unsigned int> &vertex_indices)
{
+ AssertDimension(vertex_indices.size(), GeometryInfo<1>::vertices_per_cell);
+
return all_vertices[vertex_indices[1]][0] -
all_vertices[vertex_indices[0]][0];
}
template <>
double
- cell_measure<2>(
- const std::vector<Point<2>> &all_vertices,
- const unsigned int (&vertex_indices)[GeometryInfo<2>::vertices_per_cell])
+ cell_measure<2>(const std::vector<Point<2>> & all_vertices,
+ const ArrayView<const unsigned int> &vertex_indices)
{
+ AssertDimension(vertex_indices.size(), GeometryInfo<2>::vertices_per_cell);
+
/*
Get the computation of the measure by this little Maple script. We
use the blinear mapping of the unit quad to the real quad. However,
template <>
double
- cell_measure<3>(
- const std::vector<Point<3>> &all_vertices,
- const unsigned int (&vertex_indices)[GeometryInfo<3>::vertices_per_cell])
+ cell_measure<3>(const std::vector<Point<3>> & all_vertices,
+ const ArrayView<const unsigned int> &vertex_indices)
{
+ AssertDimension(vertex_indices.size(), GeometryInfo<3>::vertices_per_cell);
// note that this is the
// cell_measure based on the new
// deal.II numbering. When called
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2018 - 2020 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.
+//
+// ---------------------------------------------------------------------
+
+
+
+// check GridTools::cell_measure()
+
+#include <deal.II/grid/grid_tools.h>
+
+#include "../tests.h"
+
+
+template <int dim>
+void
+test()
+{}
+
+template <>
+void
+test<1>()
+{
+ std::vector<Point<1>> points;
+ points.emplace_back(0.);
+ points.emplace_back(3.);
+
+ unsigned int indices1[GeometryInfo<1>::vertices_per_cell] = {0, 1};
+ std::vector<unsigned int> indices2 = {0, 1};
+
+ deallog << "1d: " << GridTools::cell_measure(points, indices1) << " "
+ << GridTools::cell_measure(points, indices2) << std::endl;
+}
+
+template <>
+void
+test<2>()
+{
+ std::vector<Point<2>> points;
+ points.emplace_back(0., 0.);
+ points.emplace_back(1., 0.);
+ points.emplace_back(0., 2.);
+ points.emplace_back(1., 2.);
+
+ unsigned int indices1[GeometryInfo<2>::vertices_per_cell] = {0, 1, 2, 3};
+ std::vector<unsigned int> indices2 = {0, 1, 2, 3};
+
+ deallog << "2d: " << GridTools::cell_measure(points, indices1) << " "
+ << GridTools::cell_measure(points, indices2) << std::endl;
+}
+
+template <>
+void
+test<3>()
+{
+ std::vector<Point<3>> points;
+ points.emplace_back(0., 0., 0.);
+ points.emplace_back(1., 0., 0.);
+ points.emplace_back(0., 2., 0.);
+ points.emplace_back(1., 2., 0.);
+ points.emplace_back(0., 0., 3.);
+ points.emplace_back(1., 0., 3.);
+ points.emplace_back(0., 2., 3.);
+ points.emplace_back(1., 2., 3.);
+
+ unsigned int indices1[GeometryInfo<3>::vertices_per_cell] = {
+ 0, 1, 2, 3, 4, 5, 6, 7};
+ std::vector<unsigned int> indices2 = {0, 1, 2, 3, 4, 5, 6, 7};
+
+ deallog << "3d: " << GridTools::cell_measure(points, indices1) << " "
+ << GridTools::cell_measure(points, indices2) << std::endl;
+}
+
+int
+main()
+{
+ initlog();
+ test<1>();
+ test<2>();
+ test<3>();
+
+ return 0;
+}