25#ifndef DUMUX_MULTIDOMAIN_EMBEDDED_LOCAL_REFINEMENT_QUADRATURE_HH
26#define DUMUX_MULTIDOMAIN_EMBEDDED_LOCAL_REFINEMENT_QUADRATURE_HH
33#include <dune/common/exceptions.hh>
52template<
class Geometry,
class IndicatorFunction>
55 using IndicatorTriple = std::array<std::size_t, 3>;
56 using GlobalPosition =
typename Geometry::GlobalCoordinate;
57 using LocalPosition =
typename Geometry::LocalCoordinate;
58 using Corners = std::array<GlobalPosition, 3>;
62 LocalPosition localPos_;
64 std::size_t indicator_;
67 QuadPoint(LocalPosition&& localPos,
double weight, std::size_t i)
68 : localPos_(weight*std::move(localPos))
73 void add(LocalPosition&& localPos,
double weight)
75 localPos_ += weight*localPos;
84 const LocalPosition& position()
const {
return localPos_; }
85 double weight()
const {
return weight_; }
86 std::size_t indicator()
const {
return indicator_; }
93 Triangle(
const Corners& c) : corners_(c) {}
94 Triangle(Corners&& c) : corners_(std::move(c)) {}
95 const GlobalPosition& operator[](std::size_t i)
const {
return corners_[i]; }
97 GlobalPosition center()
const
99 GlobalPosition center(0.0);
100 for (
const auto& c : corners_)
102 center /= corners_.size();
108 const auto ab = corners_[1] - corners_[0];
109 const auto ac = corners_[2] - corners_[0];
117 , maxLevel_(maxLevel)
120 static constexpr int dim = Geometry::mydimension;
121 static_assert(dim == 2,
"Only triangles are supported so far");
123 if (geo.corners() != (dim+1))
124 DUNE_THROW(Dune::InvalidStateException,
"Only simplex geometries are allowed");
127 const auto tri = Triangle{{ geo.corner(0), geo.corner(1), geo.corner(2) }};
128 const auto triple = IndicatorTriple{{ ind_(tri[0]), ind_(tri[1]), ind_(tri[2]) }};
129 volume_ = tri.volume();
132 addSimplex_(tri, 0, triple, triple);
134 for (
auto& qp : qps_)
138 auto begin()
const {
return qps_.begin(); }
139 auto end()
const {
return qps_.end(); }
140 auto size()
const {
return qps_.size(); }
143 void addSimplex_(
const Triangle& tri, std::size_t level,
const IndicatorTriple& triple,
const IndicatorTriple& childTriple)
146 if (std::all_of(childTriple.begin(), childTriple.end(), [a0=childTriple[0]](
auto a){ return (a == a0); }))
149 auto it = std::find_if(qps_.begin(), qps_.end(), [&](
const auto& qp){ return (qp.indicator() == childTriple[0]); });
150 if (it != qps_.end())
151 it->add(geometry_.local(tri.center()), 0.5*tri.volume()/volume_);
153 qps_.emplace_back(geometry_.local(tri.center()), 0.5*tri.volume()/volume_, childTriple[0]);
157 else if (level == maxLevel_)
159 for (
int i = 0; i < 3; ++i)
162 auto it = std::find_if(qps_.begin(), qps_.end(), [&](
const auto& qp){ return (qp.indicator() == childTriple[i]); });
163 if (it != qps_.end())
164 it->add(geometry_.local(tri[i]), 0.5*tri.volume()/volume_/3.0);
166 qps_.emplace_back(geometry_.local(tri[i]), 0.5*tri.volume()/volume_/3.0, childTriple[i]);
173 const auto children = refine(tri);
174 for (
const auto& c : children)
177 const auto grandChildTriple = IndicatorTriple{{ ind_(c[0]), ind_(c[1]), ind_(c[2]) }};
178 addSimplex_(c, level+1, triple, grandChildTriple);
183 std::array<Triangle, 4> refine(
const Triangle& tri)
const
185 const auto p01 = 0.5*(tri[0] + tri[1]);
186 const auto p12 = 0.5*(tri[1] + tri[2]);
187 const auto p20 = 0.5*(tri[2] + tri[0]);
190 {{ tri[0], p01, p20 }},
191 {{ tri[1], p01, p12 }},
192 {{ tri[2], p12, p20 }},
197 const IndicatorFunction &ind_;
198 std::size_t maxLevel_;
199 const Geometry& geometry_;
202 std::vector<QuadPoint> qps_;
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:654
Scalar volume(Shape shape, Scalar inscribedRadius)
Returns the volume of a given geometry based on the inscribed radius.
Definition: poreproperties.hh:73
Definition: localrefinementquadrature.hh:54
auto end() const
Definition: localrefinementquadrature.hh:139
auto size() const
Definition: localrefinementquadrature.hh:140
LocalRefinementSimplexQuadrature(const Geometry &geo, std::size_t maxLevel, const IndicatorFunction &i)
Definition: localrefinementquadrature.hh:115
auto begin() const
Definition: localrefinementquadrature.hh:138