From: David Wells Date: Fri, 12 Aug 2022 13:01:25 +0000 (-0400) Subject: Pass ReferenceCell by reference in QProjector. X-Git-Tag: v9.5.0-rc1~1043^2~1 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=18d870611965b8581e4f347c9906825dca79f5ec;p=dealii.git Pass ReferenceCell by reference in QProjector. This lets us use a forward declaration. --- diff --git a/include/deal.II/base/qprojector.h b/include/deal.II/base/qprojector.h index a30dd2b3aa..80d125d7cc 100644 --- a/include/deal.II/base/qprojector.h +++ b/include/deal.II/base/qprojector.h @@ -22,12 +22,13 @@ #include #include -#include - #include DEAL_II_NAMESPACE_OPEN +#ifndef DOXYGEN +class ReferenceCell; +#endif /** * @addtogroup Quadrature @@ -93,7 +94,7 @@ public: * doc for this class. */ static void - project_to_face(const ReferenceCell reference_cell, + project_to_face(const ReferenceCell & reference_cell, const SubQuadrature & quadrature, const unsigned int face_no, std::vector> &q_points); @@ -104,7 +105,7 @@ public: * the general doc for this class. */ static Quadrature - project_to_face(const ReferenceCell reference_cell, + project_to_face(const ReferenceCell &reference_cell, const SubQuadrature &quadrature, const unsigned int face_no); @@ -118,7 +119,7 @@ public: * same as those of the original rule. */ static void - project_to_subface(const ReferenceCell reference_cell, + project_to_subface(const ReferenceCell & reference_cell, const SubQuadrature & quadrature, const unsigned int face_no, const unsigned int subface_no, @@ -140,7 +141,7 @@ public: * version of this function that takes the reference cell type instead. */ static Quadrature - project_to_subface(const ReferenceCell reference_cell, + project_to_subface(const ReferenceCell & reference_cell, const SubQuadrature & quadrature, const unsigned int face_no, const unsigned int subface_no, @@ -164,7 +165,7 @@ public: * each face, in order to cope possibly different orientations of the mesh. */ static Quadrature - project_to_all_faces(const ReferenceCell reference_cell, + project_to_all_faces(const ReferenceCell & reference_cell, const hp::QCollection &quadrature); /** @@ -172,7 +173,7 @@ public: * formula on all faces. */ static Quadrature - project_to_all_faces(const ReferenceCell reference_cell, + project_to_all_faces(const ReferenceCell & reference_cell, const Quadrature &quadrature); /** @@ -189,7 +190,7 @@ public: * in FESubfaceValues. */ static Quadrature - project_to_all_subfaces(const ReferenceCell reference_cell, + project_to_all_subfaces(const ReferenceCell &reference_cell, const SubQuadrature &quadrature); /** @@ -203,7 +204,7 @@ public: * GeometryInfo::children_per_cell. */ static Quadrature - project_to_child(const ReferenceCell reference_cell, + project_to_child(const ReferenceCell & reference_cell, const Quadrature &quadrature, const unsigned int child_no); @@ -217,7 +218,7 @@ public: * refinement of the cell. */ static Quadrature - project_to_all_children(const ReferenceCell reference_cell, + project_to_all_children(const ReferenceCell & reference_cell, const Quadrature &quadrature); /** @@ -225,7 +226,7 @@ public: * connecting the points p1 and p2. */ static Quadrature - project_to_line(const ReferenceCell reference_cell, + project_to_line(const ReferenceCell &reference_cell, const Quadrature<1> &quadrature, const Point & p1, const Point & p2); @@ -270,19 +271,19 @@ public: * onto the faces) has. */ static DataSetDescriptor - face(const ReferenceCell reference_cell, - const unsigned int face_no, - const bool face_orientation, - const bool face_flip, - const bool face_rotation, - const unsigned int n_quadrature_points); + face(const ReferenceCell &reference_cell, + const unsigned int face_no, + const bool face_orientation, + const bool face_flip, + const bool face_rotation, + const unsigned int n_quadrature_points); /** * Like the above function but taking a quadrature collection, enabling * that each face might have different number of quadrature points. */ static DataSetDescriptor - face(const ReferenceCell reference_cell, + face(const ReferenceCell & reference_cell, const unsigned int face_no, const bool face_orientation, const bool face_flip, @@ -302,7 +303,7 @@ public: * Through the last argument anisotropic refinement can be respected. */ static DataSetDescriptor - subface(const ReferenceCell reference_cell, + subface(const ReferenceCell & reference_cell, const unsigned int face_no, const unsigned int subface_no, const bool face_orientation, @@ -374,7 +375,7 @@ inline QProjector::DataSetDescriptor::operator unsigned int() const template Quadrature inline QProjector::project_to_all_faces( - const ReferenceCell reference_cell, + const ReferenceCell & reference_cell, const Quadrature &quadrature) { return project_to_all_faces(reference_cell, @@ -389,32 +390,32 @@ Quadrature inline QProjector::project_to_all_faces( template <> void -QProjector<1>::project_to_face(const ReferenceCell reference_cell, +QProjector<1>::project_to_face(const ReferenceCell &reference_cell, const Quadrature<0> &, const unsigned int, std::vector> &); template <> void -QProjector<2>::project_to_face(const ReferenceCell reference_cell, +QProjector<2>::project_to_face(const ReferenceCell & reference_cell, const Quadrature<1> & quadrature, const unsigned int face_no, std::vector> &q_points); template <> void -QProjector<3>::project_to_face(const ReferenceCell reference_cell, +QProjector<3>::project_to_face(const ReferenceCell & reference_cell, const Quadrature<2> & quadrature, const unsigned int face_no, std::vector> &q_points); template <> Quadrature<1> -QProjector<1>::project_to_all_faces(const ReferenceCell reference_cell, +QProjector<1>::project_to_all_faces(const ReferenceCell & reference_cell, const hp::QCollection<0> &quadrature); template <> void -QProjector<1>::project_to_subface(const ReferenceCell reference_cell, +QProjector<1>::project_to_subface(const ReferenceCell &reference_cell, const Quadrature<0> &, const unsigned int, const unsigned int, @@ -422,7 +423,7 @@ QProjector<1>::project_to_subface(const ReferenceCell reference_cell, const RefinementCase<0> &); template <> void -QProjector<2>::project_to_subface(const ReferenceCell reference_cell, +QProjector<2>::project_to_subface(const ReferenceCell & reference_cell, const Quadrature<1> & quadrature, const unsigned int face_no, const unsigned int subface_no, @@ -430,7 +431,7 @@ QProjector<2>::project_to_subface(const ReferenceCell reference_cell, const RefinementCase<1> &); template <> void -QProjector<3>::project_to_subface(const ReferenceCell reference_cell, +QProjector<3>::project_to_subface(const ReferenceCell & reference_cell, const Quadrature<2> & quadrature, const unsigned int face_no, const unsigned int subface_no, @@ -439,7 +440,7 @@ QProjector<3>::project_to_subface(const ReferenceCell reference_cell, template <> Quadrature<1> -QProjector<1>::project_to_all_subfaces(const ReferenceCell reference_cell, +QProjector<1>::project_to_all_subfaces(const ReferenceCell &reference_cell, const Quadrature<0> &quadrature); diff --git a/source/base/qprojector.cc b/source/base/qprojector.cc index 5b5e763010..3c95d1c1cf 100644 --- a/source/base/qprojector.cc +++ b/source/base/qprojector.cc @@ -19,6 +19,8 @@ #include #include +#include + DEAL_II_NAMESPACE_OPEN @@ -81,7 +83,7 @@ namespace internal template <> void -QProjector<1>::project_to_face(const ReferenceCell reference_cell, +QProjector<1>::project_to_face(const ReferenceCell &reference_cell, const Quadrature<0> &, const unsigned int face_no, std::vector> &q_points) @@ -100,7 +102,7 @@ QProjector<1>::project_to_face(const ReferenceCell reference_cell, template <> void -QProjector<2>::project_to_face(const ReferenceCell reference_cell, +QProjector<2>::project_to_face(const ReferenceCell & reference_cell, const Quadrature<1> & quadrature, const unsigned int face_no, std::vector> &q_points) @@ -162,7 +164,7 @@ QProjector<2>::project_to_face(const ReferenceCell reference_cell, template <> void -QProjector<3>::project_to_face(const ReferenceCell reference_cell, +QProjector<3>::project_to_face(const ReferenceCell & reference_cell, const Quadrature<2> & quadrature, const unsigned int face_no, std::vector> &q_points) @@ -212,7 +214,7 @@ QProjector<3>::project_to_face(const ReferenceCell reference_cell, template <> void -QProjector<1>::project_to_subface(const ReferenceCell reference_cell, +QProjector<1>::project_to_subface(const ReferenceCell &reference_cell, const Quadrature<0> &, const unsigned int face_no, const unsigned int, @@ -233,7 +235,7 @@ QProjector<1>::project_to_subface(const ReferenceCell reference_cell, template <> void -QProjector<2>::project_to_subface(const ReferenceCell reference_cell, +QProjector<2>::project_to_subface(const ReferenceCell & reference_cell, const Quadrature<1> & quadrature, const unsigned int face_no, const unsigned int subface_no, @@ -377,7 +379,7 @@ QProjector<2>::project_to_subface(const ReferenceCell reference_cell, template <> void -QProjector<3>::project_to_subface(const ReferenceCell reference_cell, +QProjector<3>::project_to_subface(const ReferenceCell & reference_cell, const Quadrature<2> & quadrature, const unsigned int face_no, const unsigned int subface_no, @@ -467,7 +469,7 @@ QProjector<3>::project_to_subface(const ReferenceCell reference_cell, template <> Quadrature<1> -QProjector<1>::project_to_all_faces(const ReferenceCell reference_cell, +QProjector<1>::project_to_all_faces(const ReferenceCell & reference_cell, const hp::QCollection<0> &quadrature) { AssertDimension(quadrature.size(), 1); @@ -514,7 +516,7 @@ QProjector<1>::project_to_all_faces(const ReferenceCell reference_cell, template <> Quadrature<2> -QProjector<2>::project_to_all_faces(const ReferenceCell reference_cell, +QProjector<2>::project_to_all_faces(const ReferenceCell & reference_cell, const hp::QCollection<1> &quadrature) { if (reference_cell == ReferenceCells::Triangle) @@ -638,7 +640,7 @@ QProjector<2>::project_to_all_faces(const ReferenceCell reference_cell, template <> Quadrature<3> -QProjector<3>::project_to_all_faces(const ReferenceCell reference_cell, +QProjector<3>::project_to_all_faces(const ReferenceCell & reference_cell, const hp::QCollection<2> &quadrature) { const auto support_points_tri = @@ -945,7 +947,7 @@ QProjector<3>::project_to_all_faces(const ReferenceCell reference_cell, template <> Quadrature<1> -QProjector<1>::project_to_all_subfaces(const ReferenceCell reference_cell, +QProjector<1>::project_to_all_subfaces(const ReferenceCell &reference_cell, const Quadrature<0> &quadrature) { Assert(reference_cell == ReferenceCells::Line, ExcNotImplemented()); @@ -992,7 +994,7 @@ QProjector<1>::project_to_all_subfaces(const ReferenceCell reference_cell, template <> Quadrature<2> -QProjector<2>::project_to_all_subfaces(const ReferenceCell reference_cell, +QProjector<2>::project_to_all_subfaces(const ReferenceCell &reference_cell, const SubQuadrature &quadrature) { if (reference_cell == ReferenceCells::Triangle || @@ -1043,7 +1045,7 @@ QProjector<2>::project_to_all_subfaces(const ReferenceCell reference_cell, template <> Quadrature<3> -QProjector<3>::project_to_all_subfaces(const ReferenceCell reference_cell, +QProjector<3>::project_to_all_subfaces(const ReferenceCell &reference_cell, const SubQuadrature &quadrature) { if (reference_cell == ReferenceCells::Triangle || @@ -1127,7 +1129,7 @@ QProjector<3>::project_to_all_subfaces(const ReferenceCell reference_cell, template Quadrature -QProjector::project_to_child(const ReferenceCell reference_cell, +QProjector::project_to_child(const ReferenceCell & reference_cell, const Quadrature &quadrature, const unsigned int child_no) { @@ -1159,7 +1161,7 @@ QProjector::project_to_child(const ReferenceCell reference_cell, template Quadrature -QProjector::project_to_all_children(const ReferenceCell reference_cell, +QProjector::project_to_all_children(const ReferenceCell & reference_cell, const Quadrature &quadrature) { Assert(reference_cell == ReferenceCells::get_hypercube(), @@ -1191,7 +1193,7 @@ QProjector::project_to_all_children(const ReferenceCell reference_cell, template Quadrature -QProjector::project_to_line(const ReferenceCell reference_cell, +QProjector::project_to_line(const ReferenceCell &reference_cell, const Quadrature<1> &quadrature, const Point & p1, const Point & p2) @@ -1219,11 +1221,11 @@ QProjector::project_to_line(const ReferenceCell reference_cell, template typename QProjector::DataSetDescriptor -QProjector::DataSetDescriptor::face(const ReferenceCell reference_cell, - const unsigned int face_no, - const bool face_orientation, - const bool face_flip, - const bool face_rotation, +QProjector::DataSetDescriptor::face(const ReferenceCell &reference_cell, + const unsigned int face_no, + const bool face_orientation, + const bool face_flip, + const bool face_rotation, const unsigned int n_quadrature_points) { if (reference_cell == ReferenceCells::Triangle || @@ -1309,7 +1311,7 @@ QProjector::DataSetDescriptor::face(const ReferenceCell reference_cell, template typename QProjector::DataSetDescriptor QProjector::DataSetDescriptor::face( - const ReferenceCell reference_cell, + const ReferenceCell & reference_cell, const unsigned int face_no, const bool face_orientation, const bool face_flip, @@ -1446,9 +1448,9 @@ QProjector::DataSetDescriptor::face( template <> QProjector<1>::DataSetDescriptor QProjector<1>::DataSetDescriptor::subface( - const ReferenceCell reference_cell, - const unsigned int face_no, - const unsigned int subface_no, + const ReferenceCell &reference_cell, + const unsigned int face_no, + const unsigned int subface_no, const bool, const bool, const bool, @@ -1471,9 +1473,9 @@ QProjector<1>::DataSetDescriptor::subface( template <> QProjector<2>::DataSetDescriptor QProjector<2>::DataSetDescriptor::subface( - const ReferenceCell reference_cell, - const unsigned int face_no, - const unsigned int subface_no, + const ReferenceCell &reference_cell, + const unsigned int face_no, + const unsigned int subface_no, const bool, const bool, const bool, @@ -1496,7 +1498,7 @@ QProjector<2>::DataSetDescriptor::subface( template <> QProjector<3>::DataSetDescriptor QProjector<3>::DataSetDescriptor::subface( - const ReferenceCell reference_cell, + const ReferenceCell & reference_cell, const unsigned int face_no, const unsigned int subface_no, const bool face_orientation, @@ -1753,7 +1755,7 @@ QProjector<3>::DataSetDescriptor::subface( template Quadrature -QProjector::project_to_face(const ReferenceCell reference_cell, +QProjector::project_to_face(const ReferenceCell &reference_cell, const SubQuadrature &quadrature, const unsigned int face_no) { @@ -1770,7 +1772,7 @@ QProjector::project_to_face(const ReferenceCell reference_cell, template Quadrature -QProjector::project_to_subface(const ReferenceCell reference_cell, +QProjector::project_to_subface(const ReferenceCell &reference_cell, const SubQuadrature &quadrature, const unsigned int face_no, const unsigned int subface_no,