26#ifndef DUMUX_COMMON_REFINEMENT_QUADRATURERULE_HH
27#define DUMUX_COMMON_REFINEMENT_QUADRATURERULE_HH
29#include <dune/geometry/quadraturerules.hh>
30#include <dune/geometry/virtualrefinement.hh>
39template<
typename ct,
int mydim>
44 static constexpr int dim = mydim;
54 auto& refinement = Dune::buildRefinement<dim, ct>(type, type);
55 const auto tag = Dune::refinementLevels(levels);
58 const auto numVertices = refinement.nVertices(tag);
59 std::vector<Dune::FieldVector<ct, dim>> points; points.reserve(numVertices);
60 auto vIt = refinement.vBegin(tag);
61 const auto vItEnd = refinement.vEnd(tag);
62 for (; vIt != vItEnd; ++vIt)
63 points.emplace_back(vIt.coords());
66 const auto numElements = refinement.nElements(tag);
67 this->reserve(numElements);
68 auto eIt = refinement.eBegin(tag);
69 const auto eItEnd = refinement.eEnd(tag);
70 for (; eIt != eItEnd; ++eIt)
73 const auto weight = computeVolume_(type, points, eIt.vertexIndices());
74 this->emplace_back(eIt.coords(), weight);
79 template<
class VertexIndices>
80 ct computeVolume_(Dune::GeometryType t,
const std::vector<Dune::FieldVector<ct, dim>>& points,
const VertexIndices& indices)
const
82 if (t != Dune::GeometryTypes::simplex(2))
83 DUNE_THROW(Dune::NotImplemented,
"Only implemented for 2d simplices");
85 const auto ab = points[indices[1]] - points[indices[0]];
86 const auto ac = points[indices[2]] - points[indices[0]];
Define some often used mathematical functions.
Dune::FieldVector< Scalar, 3 > crossProduct(const Dune::FieldVector< Scalar, 3 > &vec1, const Dune::FieldVector< Scalar, 3 > &vec2)
Cross product of two vectors in three-dimensional Euclidean space.
Definition: math.hh:631
Definition: common/pdesolver.hh:35
A "quadrature" based on virtual refinement.
Definition: refinementquadraturerule.hh:41
~RefinementQuadratureRule() final
Definition: refinementquadraturerule.hh:49
RefinementQuadratureRule(Dune::GeometryType type, int levels)
Definition: refinementquadraturerule.hh:51
static constexpr int dim
The space dimension.
Definition: refinementquadraturerule.hh:44
static constexpr int highest_order
The highest quadrature order available.
Definition: refinementquadraturerule.hh:47