diff --git a/src/analysis/QuadratureManager.hpp b/src/analysis/QuadratureManager.hpp index b0f81c543c0612181cc5728b2c7ef36fc1e25515..04655f291f8ff75e8c1207c4b6171a7b7f461a05 100644 --- a/src/analysis/QuadratureManager.hpp +++ b/src/analysis/QuadratureManager.hpp @@ -47,6 +47,10 @@ class QuadratureManager const QuadratureFormula<1>& getLineFormula(const IQuadratureDescriptor& quadrature_descriptor) const { + if (quadrature_descriptor.degree() > maxLineDegree(quadrature_descriptor.type())) { + throw NormalError(::name(quadrature_descriptor.type()) + " quadrature formulae handle degrees up to " + + std::to_string(maxLineDegree(quadrature_descriptor.type())) + " on lines"); + } switch (quadrature_descriptor.type()) { case QuadratureType::Gauss: case QuadratureType::GaussLegendre: { @@ -56,28 +60,40 @@ class QuadratureManager case QuadratureType::GaussLobatto: { return m_line_gauss_lobatto_formula_list[quadrature_descriptor.degree() / 2]; } + // LCOV_EXCL_START default: { - throw UnexpectedError("invalid quadrature type"); + throw UnexpectedError(::name(quadrature_descriptor.type()) + " is not defined on lines"); } + // LCOV_EXCL_STOP } } const QuadratureFormula<2>& getTriangleFormula(const IQuadratureDescriptor& quadrature_descriptor) const { + if (quadrature_descriptor.degree() > maxTriangleDegree(quadrature_descriptor.type())) { + throw NormalError(::name(quadrature_descriptor.type()) + " quadrature formulae handle degrees up to " + + std::to_string(maxTriangleDegree(quadrature_descriptor.type())) + " on triangles"); + } switch (quadrature_descriptor.type()) { case QuadratureType::Gauss: { return m_triangle_gauss_formula_list[quadrature_descriptor.degree() - 1]; } + // LCOV_EXCL_START default: { - throw UnexpectedError(quadrature_descriptor.name() + " is not defined on triangles"); + throw UnexpectedError(::name(quadrature_descriptor.type()) + " is not defined on triangles"); } + // LCOV_EXCL_STOP } } const QuadratureFormula<2>& getSquareFormula(const IQuadratureDescriptor& quadrature_descriptor) const { + if (quadrature_descriptor.degree() > maxSquareDegree(quadrature_descriptor.type())) { + throw NormalError(::name(quadrature_descriptor.type()) + " quadrature formulae handle degrees up to " + + std::to_string(maxSquareDegree(quadrature_descriptor.type())) + " on squares"); + } switch (quadrature_descriptor.type()) { case QuadratureType::Gauss: { return m_square_gauss_formula_list[quadrature_descriptor.degree() / 2]; @@ -88,54 +104,78 @@ class QuadratureManager case QuadratureType::GaussLobatto: { return m_square_gauss_lobatto_formula_list[quadrature_descriptor.degree() / 2]; } + // LCOV_EXCL_START default: { - throw UnexpectedError("invalid quadrature type"); + throw UnexpectedError(::name(quadrature_descriptor.type()) + " is not defined on squares"); } + // LCOV_EXCL_STOP } } const QuadratureFormula<3>& getTetrahedronFormula(const IQuadratureDescriptor& quadrature_descriptor) const { + if (quadrature_descriptor.degree() > maxTetrahedronDegree(quadrature_descriptor.type())) { + throw NormalError(::name(quadrature_descriptor.type()) + " quadrature formulae handle degrees up to " + + std::to_string(maxTetrahedronDegree(quadrature_descriptor.type())) + " on tetrahedra"); + } switch (quadrature_descriptor.type()) { case QuadratureType::Gauss: { return m_tetrahedron_gauss_formula_list[quadrature_descriptor.degree() - 1]; } + // LCOV_EXCL_START default: { - throw UnexpectedError(quadrature_descriptor.name() + " is not defined on tetrahedron"); + throw UnexpectedError(::name(quadrature_descriptor.type()) + " is not defined on tetrahedra"); } + // LCOV_EXCL_STOP } } const QuadratureFormula<3>& getPrismFormula(const IQuadratureDescriptor& quadrature_descriptor) const { + if (quadrature_descriptor.degree() > maxPrismDegree(quadrature_descriptor.type())) { + throw NormalError(::name(quadrature_descriptor.type()) + " quadrature formulae handle degrees up to " + + std::to_string(maxPrismDegree(quadrature_descriptor.type())) + " on prisms"); + } switch (quadrature_descriptor.type()) { case QuadratureType::Gauss: { return m_prism_gauss_formula_list[quadrature_descriptor.degree() - 1]; } + // LCOV_EXCL_START default: { - throw UnexpectedError(quadrature_descriptor.name() + " is not defined on prism"); + throw UnexpectedError(::name(quadrature_descriptor.type()) + " is not defined on prisms"); } + // LCOV_EXCL_STOP } } const QuadratureFormula<3>& getPyramidFormula(const IQuadratureDescriptor& quadrature_descriptor) const { + if (quadrature_descriptor.degree() > maxPyramidDegree(quadrature_descriptor.type())) { + throw NormalError(::name(quadrature_descriptor.type()) + " quadrature formulae handle degrees up to " + + std::to_string(maxPyramidDegree(quadrature_descriptor.type())) + " on pyramids"); + } switch (quadrature_descriptor.type()) { case QuadratureType::Gauss: { return m_pyramid_gauss_formula_list[quadrature_descriptor.degree() - 1]; } + // LCOV_EXCL_START default: { - throw UnexpectedError(quadrature_descriptor.name() + " is not defined on pyramid"); + throw UnexpectedError(::name(quadrature_descriptor.type()) + " is not defined on pyramid"); } + // LCOV_EXCL_STOP } } const QuadratureFormula<3>& getCubeFormula(const IQuadratureDescriptor& quadrature_descriptor) const { + if (quadrature_descriptor.degree() > maxCubeDegree(quadrature_descriptor.type())) { + throw NormalError(::name(quadrature_descriptor.type()) + " quadrature formulae handle degrees up to " + + std::to_string(maxCubeDegree(quadrature_descriptor.type())) + " on cubes"); + } switch (quadrature_descriptor.type()) { case QuadratureType::Gauss: { return m_cube_gauss_formula_list[quadrature_descriptor.degree() / 2]; @@ -147,9 +187,11 @@ class QuadratureManager case QuadratureType::GaussLobatto: { return m_cube_gauss_lobatto_formula_list[quadrature_descriptor.degree() / 2]; } + // LCOV_EXCL_START default: { - throw UnexpectedError("invalid quadrature type"); + throw UnexpectedError(::name(quadrature_descriptor.type()) + " is not defined on cubes"); } + // LCOV_EXCL_STOP } } @@ -164,9 +206,11 @@ class QuadratureManager case QuadratureType::GaussLobatto: { return m_line_gauss_lobatto_formula_list.size() * 2 - 1; } + // LCOV_EXCL_START default: { throw UnexpectedError("invalid quadrature type"); } + // LCOV_EXCL_STOP } } @@ -183,9 +227,11 @@ class QuadratureManager case QuadratureType::GaussLobatto: { return m_square_gauss_lobatto_formula_list.size() * 2 - 1; } + // LCOV_EXCL_START default: { throw UnexpectedError("invalid quadrature type"); } + // LCOV_EXCL_STOP } } @@ -196,9 +242,11 @@ class QuadratureManager case QuadratureType::Gauss: { return m_triangle_gauss_formula_list.size(); } + // LCOV_EXCL_START default: { - throw UnexpectedError(name(type) + " is not defined on triangle"); + throw UnexpectedError(::name(type) + " is not defined on triangle"); } + // LCOV_EXCL_STOP } } @@ -215,9 +263,11 @@ class QuadratureManager case QuadratureType::GaussLobatto: { return m_cube_gauss_lobatto_formula_list.size() * 2 - 1; } + // LCOV_EXCL_START default: { throw UnexpectedError("invalid quadrature type"); } + // LCOV_EXCL_STOP } } @@ -228,9 +278,11 @@ class QuadratureManager case QuadratureType::Gauss: { return m_prism_gauss_formula_list.size(); } + // LCOV_EXCL_START default: { throw UnexpectedError(::name(type) + " is not defined on prism"); } + // LCOV_EXCL_STOP } } @@ -241,9 +293,11 @@ class QuadratureManager case QuadratureType::Gauss: { return m_pyramid_gauss_formula_list.size(); } + // LCOV_EXCL_START default: { throw UnexpectedError(::name(type) + " is not defined on pyramid"); } + // LCOV_EXCL_STOP } } @@ -254,9 +308,11 @@ class QuadratureManager case QuadratureType::Gauss: { return m_tetrahedron_gauss_formula_list.size(); } + // LCOV_EXCL_START default: { throw UnexpectedError(::name(type) + " is not defined on tetrahedron"); } + // LCOV_EXCL_STOP } } diff --git a/tests/test_CubeGaussQuadrature.cpp b/tests/test_CubeGaussQuadrature.cpp index 1db30eb178be24be223f63b5f6a68aa441261909..d4221fe322dc6dc4f898ac2563dc14d463b2e71f 100644 --- a/tests/test_CubeGaussQuadrature.cpp +++ b/tests/test_CubeGaussQuadrature.cpp @@ -488,6 +488,10 @@ TEST_CASE("CubeGaussQuadrature", "[analysis]") SECTION("max implemented degree") { REQUIRE(QuadratureManager::instance().maxCubeDegree(QuadratureType::Gauss) == CubeGaussQuadrature::max_degree); + REQUIRE_THROWS_WITH(QuadratureManager::instance().getCubeFormula( + GaussQuadratureDescriptor(CubeGaussQuadrature ::max_degree + 1)), + "error: Gauss quadrature formulae handle degrees up to " + + std::to_string(CubeGaussQuadrature ::max_degree) + " on cubes"); } SECTION("Access functions") diff --git a/tests/test_PrismGaussQuadrature.cpp b/tests/test_PrismGaussQuadrature.cpp index ceff287e5ac5c823c9f427ba6801babfeaa47481..000edeed6a40556d702bd440a26aab90f37b7e5d 100644 --- a/tests/test_PrismGaussQuadrature.cpp +++ b/tests/test_PrismGaussQuadrature.cpp @@ -682,5 +682,9 @@ TEST_CASE("PrismGaussQuadrature", "[analysis]") SECTION("max implemented degree") { REQUIRE(QuadratureManager::instance().maxPrismDegree(QuadratureType::Gauss) == PrismGaussQuadrature::max_degree); + REQUIRE_THROWS_WITH(QuadratureManager::instance().getPrismFormula( + GaussQuadratureDescriptor(PrismGaussQuadrature ::max_degree + 1)), + "error: Gauss quadrature formulae handle degrees up to " + + std::to_string(PrismGaussQuadrature ::max_degree) + " on prisms"); } } diff --git a/tests/test_PyramidGaussQuadrature.cpp b/tests/test_PyramidGaussQuadrature.cpp index e648c96fb13df144945e7c73cdf381545ff1169b..ea022b1360908ccd3f270d3afda53b710a486488 100644 --- a/tests/test_PyramidGaussQuadrature.cpp +++ b/tests/test_PyramidGaussQuadrature.cpp @@ -566,5 +566,9 @@ TEST_CASE("PyramidGaussQuadrature", "[analysis]") { REQUIRE(QuadratureManager::instance().maxPyramidDegree(QuadratureType::Gauss) == PyramidGaussQuadrature::max_degree); + REQUIRE_THROWS_WITH(QuadratureManager::instance().getPyramidFormula( + GaussQuadratureDescriptor(PyramidGaussQuadrature ::max_degree + 1)), + "error: Gauss quadrature formulae handle degrees up to " + + std::to_string(PyramidGaussQuadrature ::max_degree) + " on pyramids"); } } diff --git a/tests/test_SquareGaussQuadrature.cpp b/tests/test_SquareGaussQuadrature.cpp index 5e6b9d6e5ee066cd8db48c6a2297bfc3d12cfcd8..74f442a3154e23e4436c9b170a9649e72a61d4ae 100644 --- a/tests/test_SquareGaussQuadrature.cpp +++ b/tests/test_SquareGaussQuadrature.cpp @@ -459,6 +459,10 @@ TEST_CASE("SquareGaussQuadrature", "[analysis]") SECTION("max implemented degree") { REQUIRE(QuadratureManager::instance().maxSquareDegree(QuadratureType::Gauss) == SquareGaussQuadrature::max_degree); + REQUIRE_THROWS_WITH(QuadratureManager::instance().getSquareFormula( + GaussQuadratureDescriptor(SquareGaussQuadrature ::max_degree + 1)), + "error: Gauss quadrature formulae handle degrees up to " + + std::to_string(SquareGaussQuadrature ::max_degree) + " on squares"); } SECTION("Access functions") diff --git a/tests/test_TensorialGaussLegendreQuadrature.cpp b/tests/test_TensorialGaussLegendreQuadrature.cpp index aa8708ec610e0c43fa9ad2ded1f50c78e5461b91..342a4fbaa9ef8c54e274f7ad65800aee30986559 100644 --- a/tests/test_TensorialGaussLegendreQuadrature.cpp +++ b/tests/test_TensorialGaussLegendreQuadrature.cpp @@ -3,6 +3,7 @@ #include <catch2/matchers/catch_matchers_all.hpp> #include <analysis/GaussLegendreQuadratureDescriptor.hpp> +#include <analysis/GaussQuadratureDescriptor.hpp> #include <analysis/QuadratureManager.hpp> #include <analysis/TensorialGaussLegendreQuadrature.hpp> #include <utils/Exceptions.hpp> @@ -519,6 +520,15 @@ TEST_CASE("TensorialGaussLegendreQuadrature", "[analysis]") { REQUIRE(QuadratureManager::instance().maxLineDegree(QuadratureType::GaussLegendre) == TensorialGaussLegendreQuadrature<1>::max_degree); + REQUIRE_THROWS_WITH(QuadratureManager::instance().getLineFormula( + GaussLegendreQuadratureDescriptor(TensorialGaussLegendreQuadrature<1>::max_degree + 1)), + "error: Gauss-Legendre quadrature formulae handle degrees up to " + + std::to_string(TensorialGaussLegendreQuadrature<1>::max_degree) + " on lines"); + + REQUIRE_THROWS_WITH(QuadratureManager::instance().getLineFormula( + GaussQuadratureDescriptor(TensorialGaussLegendreQuadrature<1>::max_degree + 1)), + "error: Gauss quadrature formulae handle degrees up to " + + std::to_string(TensorialGaussLegendreQuadrature<1>::max_degree) + " on lines"); } } @@ -569,6 +579,16 @@ TEST_CASE("TensorialGaussLegendreQuadrature", "[analysis]") REQUIRE(integrate(px7y6, l7, 0, 1, 0.2, 0.8) == Catch::Approx(std::pow(0.8, 7) - std::pow(0.2, 7))); } + + SECTION("max implemented degree") + { + REQUIRE(QuadratureManager::instance().maxSquareDegree(QuadratureType::GaussLegendre) == + TensorialGaussLegendreQuadrature<2>::max_degree); + REQUIRE_THROWS_WITH(QuadratureManager::instance().getSquareFormula( + GaussLegendreQuadratureDescriptor(TensorialGaussLegendreQuadrature<2>::max_degree + 1)), + "error: Gauss-Legendre quadrature formulae handle degrees up to " + + std::to_string(TensorialGaussLegendreQuadrature<2>::max_degree) + " on squares"); + } } SECTION("3D") @@ -624,5 +644,15 @@ TEST_CASE("TensorialGaussLegendreQuadrature", "[analysis]") Catch::Approx((std::pow(0.8, 7) - std::pow(0.2, 7)) * (0.7 - -0.1) + (0.8 - 0.2) * (std::pow(0.7, 8) - std::pow(-0.1, 8)))); } + + SECTION("max implemented degree") + { + REQUIRE(QuadratureManager::instance().maxCubeDegree(QuadratureType::GaussLegendre) == + TensorialGaussLegendreQuadrature<3>::max_degree); + REQUIRE_THROWS_WITH(QuadratureManager::instance().getCubeFormula( + GaussLegendreQuadratureDescriptor(TensorialGaussLegendreQuadrature<3>::max_degree + 1)), + "error: Gauss-Legendre quadrature formulae handle degrees up to " + + std::to_string(TensorialGaussLegendreQuadrature<3>::max_degree) + " on cubes"); + } } } diff --git a/tests/test_TensorialGaussLobattoQuadrature.cpp b/tests/test_TensorialGaussLobattoQuadrature.cpp index e038b692ea6dce069e037640307e9d57864fffef..9ee4028cf0aa9b58d9f7478e91f2146f7f77a675 100644 --- a/tests/test_TensorialGaussLobattoQuadrature.cpp +++ b/tests/test_TensorialGaussLobattoQuadrature.cpp @@ -294,6 +294,10 @@ TEST_CASE("TensorialGaussLobattoQuadrature", "[analysis]") { REQUIRE(QuadratureManager::instance().maxLineDegree(QuadratureType::GaussLobatto) == TensorialGaussLobattoQuadrature<1>::max_degree); + REQUIRE_THROWS_WITH(QuadratureManager::instance().getLineFormula( + GaussLobattoQuadratureDescriptor(TensorialGaussLobattoQuadrature<1>::max_degree + 1)), + "error: Gauss-Lobatto quadrature formulae handle degrees up to " + + std::to_string(TensorialGaussLobattoQuadrature<1>::max_degree) + " on lines"); } } @@ -344,6 +348,16 @@ TEST_CASE("TensorialGaussLobattoQuadrature", "[analysis]") REQUIRE(integrate(px7y6, l7, 0, 1, 0.2, 0.8) == Catch::Approx(std::pow(0.8, 7) - std::pow(0.2, 7))); } + + SECTION("max implemented degree") + { + REQUIRE(QuadratureManager::instance().maxSquareDegree(QuadratureType::GaussLobatto) == + TensorialGaussLobattoQuadrature<2>::max_degree); + REQUIRE_THROWS_WITH(QuadratureManager::instance().getSquareFormula( + GaussLobattoQuadratureDescriptor(TensorialGaussLobattoQuadrature<2>::max_degree + 1)), + "error: Gauss-Lobatto quadrature formulae handle degrees up to " + + std::to_string(TensorialGaussLobattoQuadrature<2>::max_degree) + " on squares"); + } } SECTION("3D") @@ -399,6 +413,16 @@ TEST_CASE("TensorialGaussLobattoQuadrature", "[analysis]") Catch::Approx((std::pow(0.8, 7) - std::pow(0.2, 7)) * (0.7 - -0.1) + (0.8 - 0.2) * (std::pow(0.7, 8) - std::pow(-0.1, 8)))); } + + SECTION("max implemented degree") + { + REQUIRE(QuadratureManager::instance().maxCubeDegree(QuadratureType::GaussLobatto) == + TensorialGaussLobattoQuadrature<3>::max_degree); + REQUIRE_THROWS_WITH(QuadratureManager::instance().getCubeFormula( + GaussLobattoQuadratureDescriptor(TensorialGaussLobattoQuadrature<3>::max_degree + 1)), + "error: Gauss-Lobatto quadrature formulae handle degrees up to " + + std::to_string(TensorialGaussLobattoQuadrature<3>::max_degree) + " on cubes"); + } } SECTION("Access functions") diff --git a/tests/test_TetrahedronGaussQuadrature.cpp b/tests/test_TetrahedronGaussQuadrature.cpp index 605d1dc1fa7476719fe16ac595ff55491f5bd8e8..a412c13688596e37916fd50a522b00a6aaae8485 100644 --- a/tests/test_TetrahedronGaussQuadrature.cpp +++ b/tests/test_TetrahedronGaussQuadrature.cpp @@ -679,5 +679,9 @@ TEST_CASE("TetrahedronGaussQuadrature", "[analysis]") { REQUIRE(QuadratureManager::instance().maxTetrahedronDegree(QuadratureType::Gauss) == TetrahedronGaussQuadrature::max_degree); + REQUIRE_THROWS_WITH(QuadratureManager::instance().getTetrahedronFormula( + GaussQuadratureDescriptor(TetrahedronGaussQuadrature ::max_degree + 1)), + "error: Gauss quadrature formulae handle degrees up to " + + std::to_string(TetrahedronGaussQuadrature ::max_degree) + " on tetrahedra"); } } diff --git a/tests/test_TriangleGaussQuadrature.cpp b/tests/test_TriangleGaussQuadrature.cpp index 3d7b5dffaf2f0290d5b15dcb6796188512f62cf8..4ed18d60284861241861c45dd744a8e1b38c794b 100644 --- a/tests/test_TriangleGaussQuadrature.cpp +++ b/tests/test_TriangleGaussQuadrature.cpp @@ -624,5 +624,9 @@ TEST_CASE("TriangleGaussQuadrature", "[analysis]") { REQUIRE(QuadratureManager::instance().maxTriangleDegree(QuadratureType::Gauss) == TriangleGaussQuadrature::max_degree); + REQUIRE_THROWS_WITH(QuadratureManager::instance().getTriangleFormula( + GaussQuadratureDescriptor(TriangleGaussQuadrature::max_degree + 1)), + "error: Gauss quadrature formulae handle degrees up to " + + std::to_string(TriangleGaussQuadrature::max_degree) + " on triangles"); } }