25#ifndef DUMUX_MULTIDOMAIN_EMBEDDED_LOCAL_REFINEMENT_QUADRATURE_HH
26#define DUMUX_MULTIDOMAIN_EMBEDDED_LOCAL_REFINEMENT_QUADRATURE_HH
33#include <dune/common/exceptions.hh>
51template<
class Geometry,
class IndicatorFunction>
54 using IndicatorTriple = std::array<std::size_t, 3>;
55 using GlobalPosition =
typename Geometry::GlobalCoordinate;
56 using LocalPosition =
typename Geometry::LocalCoordinate;
57 using Corners = std::array<GlobalPosition, 3>;
61 LocalPosition localPos_;
63 std::size_t indicator_;
66 QuadPoint(LocalPosition&& localPos,
double weight, std::size_t i)
67 : localPos_(weight*std::move(localPos))
72 void add(LocalPosition&& localPos,
double weight)
74 localPos_ += weight*localPos;
83 const LocalPosition& position()
const {
return localPos_; }
84 double weight()
const {
return weight_; }
85 std::size_t indicator()
const {
return indicator_; }
92 Triangle(
const Corners& c) : corners_(c) {}
93 Triangle(Corners&& c) : corners_(std::move(c)) {}
94 const GlobalPosition& operator[](std::size_t i)
const {
return corners_[i]; }
96 GlobalPosition
center()
const
98 GlobalPosition
center(0.0);
99 for (
const auto& c : corners_)
101 center /= corners_.size();
107 const auto ab = corners_[1] - corners_[0];
108 const auto ac = corners_[2] - corners_[0];
116 , maxLevel_(maxLevel)
119 static constexpr int dim = Geometry::mydimension;
120 static_assert(dim == 2,
"Only triangles are supported so far");
122 if (geo.corners() != (dim+1))
123 DUNE_THROW(Dune::InvalidStateException,
"Only simplex geometries are allowed");
126 const auto tri = Triangle{{ geo.corner(0), geo.corner(1), geo.corner(2) }};
127 const auto triple = IndicatorTriple{{ ind_(tri[0]), ind_(tri[1]), ind_(tri[2]) }};
128 volume_ = tri.volume();
131 addSimplex_(tri, 0, triple, triple);
133 for (
auto& qp : qps_)
137 auto begin()
const {
return qps_.begin(); }
138 auto end()
const {
return qps_.end(); }
139 auto size()
const {
return qps_.size(); }
142 void addSimplex_(
const Triangle& tri, std::size_t level,
const IndicatorTriple& triple,
const IndicatorTriple& childTriple)
145 if (std::all_of(childTriple.begin(), childTriple.end(), [a0=childTriple[0]](
auto a){ return (a == a0); }))
148 auto it = std::find_if(qps_.begin(), qps_.end(), [&](
const auto& qp){ return (qp.indicator() == childTriple[0]); });
149 if (it != qps_.end())
150 it->add(geometry_.local(tri.center()), 0.5*tri.volume()/volume_);
152 qps_.emplace_back(geometry_.local(tri.center()), 0.5*tri.volume()/volume_, childTriple[0]);
156 else if (level == maxLevel_)
158 for (
int i = 0; i < 3; ++i)
161 auto it = std::find_if(qps_.begin(), qps_.end(), [&](
const auto& qp){ return (qp.indicator() == childTriple[i]); });
162 if (it != qps_.end())
163 it->add(geometry_.local(tri[i]), 0.5*tri.volume()/volume_/3.0);
165 qps_.emplace_back(geometry_.local(tri[i]), 0.5*tri.volume()/volume_/3.0, childTriple[i]);
172 const auto children = refine(tri);
173 for (
const auto& c : children)
176 const auto grandChildTriple = IndicatorTriple{{ ind_(c[0]), ind_(c[1]), ind_(c[2]) }};
177 addSimplex_(c, level+1, triple, grandChildTriple);
182 std::array<Triangle, 4> refine(
const Triangle& tri)
const
184 const auto p01 = 0.5*(tri[0] + tri[1]);
185 const auto p12 = 0.5*(tri[1] + tri[2]);
186 const auto p20 = 0.5*(tri[2] + tri[0]);
189 {{ tri[0], p01, p20 }},
190 {{ tri[1], p01, p12 }},
191 {{ tri[2], p12, p20 }},
196 const IndicatorFunction &ind_;
197 std::size_t maxLevel_;
198 const Geometry& geometry_;
201 std::vector<QuadPoint> qps_;
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:171
Corners::value_type center(const Corners &corners)
The center of a given list of corners.
Definition: center.hh:36
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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
A quadrature rule using local refinement to approximate partitioned elements.
Definition: localrefinementquadrature.hh:53
auto end() const
Definition: localrefinementquadrature.hh:138
auto size() const
Definition: localrefinementquadrature.hh:139
LocalRefinementSimplexQuadrature(const Geometry &geo, std::size_t maxLevel, const IndicatorFunction &i)
Definition: localrefinementquadrature.hh:114
auto begin() const
Definition: localrefinementquadrature.hh:137