3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
localrefinementquadrature.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_MULTIDOMAIN_EMBEDDED_LOCAL_REFINEMENT_QUADRATURE_HH
26#define DUMUX_MULTIDOMAIN_EMBEDDED_LOCAL_REFINEMENT_QUADRATURE_HH
27
28#include <array>
29#include <vector>
30#include <utility>
31#include <algorithm>
32
33#include <dune/common/exceptions.hh>
34
35namespace Dumux {
36
52template<class Geometry, class IndicatorFunction>
54{
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>;
59
60 class QuadPoint
61 {
62 LocalPosition localPos_;
63 double weight_;
64 std::size_t indicator_;
65
66 public:
67 QuadPoint(LocalPosition&& localPos, double weight, std::size_t i)
68 : localPos_(weight*std::move(localPos))
69 , weight_(weight)
70 , indicator_(i)
71 {}
72
73 void add(LocalPosition&& localPos, double weight)
74 {
75 localPos_ += weight*localPos;
76 weight_ += weight;
77 }
78
79 void finalize()
80 {
81 localPos_ /= weight_;
82 }
83
84 const LocalPosition& position() const { return localPos_; }
85 double weight() const { return weight_; }
86 std::size_t indicator() const { return indicator_; }
87 };
88
89 class Triangle
90 {
91 Corners corners_;
92 public:
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]; }
96
97 GlobalPosition center() const
98 {
99 GlobalPosition center(0.0);
100 for (const auto& c : corners_)
101 center += c;
102 center /= corners_.size();
103 return center;
104 }
105
106 auto volume() const
107 {
108 const auto ab = corners_[1] - corners_[0];
109 const auto ac = corners_[2] - corners_[0];
110 return crossProduct(ab, ac).two_norm()*0.5;
111 }
112 };
113
114public:
115 LocalRefinementSimplexQuadrature(const Geometry& geo, std::size_t maxLevel, const IndicatorFunction& i)
116 : ind_(i)
117 , maxLevel_(maxLevel)
118 , geometry_(geo)
119 {
120 static constexpr int dim = Geometry::mydimension;
121 static_assert(dim == 2, "Only triangles are supported so far");
122
123 if (geo.corners() != (dim+1))
124 DUNE_THROW(Dune::InvalidStateException, "Only simplex geometries are allowed");
125
126 // evaluate indicator
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();
130
131 // add new sub-qps recursively by adaptive refinement
132 addSimplex_(tri, 0, triple, triple);
133
134 for (auto& qp : qps_)
135 qp.finalize();
136 }
137
138 auto begin() const { return qps_.begin(); }
139 auto end() const { return qps_.end(); }
140 auto size() const { return qps_.size(); }
141
142private:
143 void addSimplex_(const Triangle& tri, std::size_t level, const IndicatorTriple& triple, const IndicatorTriple& childTriple)
144 {
145 // if all indicator are the same just add the triangle
146 if (std::all_of(childTriple.begin(), childTriple.end(), [a0=childTriple[0]](auto a){ return (a == a0); }))
147 {
148 // the volumes need to sum up to 0.5 (triangle reference element volume)
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_);
152 else
153 qps_.emplace_back(geometry_.local(tri.center()), 0.5*tri.volume()/volume_, childTriple[0]);
154 }
155
156 // if they are different but the maximum level is reached, split volume in three
157 else if (level == maxLevel_)
158 {
159 for (int i = 0; i < 3; ++i)
160 {
161 // the volumes need to sum up to 0.5 (triangle reference element volume)
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);
165 else
166 qps_.emplace_back(geometry_.local(tri[i]), 0.5*tri.volume()/volume_/3.0, childTriple[i]);
167 }
168 }
169
170 // otherwise refine and go recursively over all children
171 else
172 {
173 const auto children = refine(tri);
174 for (const auto& c : children)
175 {
176 // make sure to recompute indicator every time since new indices might come up
177 const auto grandChildTriple = IndicatorTriple{{ ind_(c[0]), ind_(c[1]), ind_(c[2]) }};
178 addSimplex_(c, level+1, triple, grandChildTriple);
179 }
180 }
181 }
182
183 std::array<Triangle, 4> refine(const Triangle& tri) const
184 {
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]);
188
189 return {{
190 {{ tri[0], p01, p20 }},
191 {{ tri[1], p01, p12 }},
192 {{ tri[2], p12, p20 }},
193 {{ p12, p01, p20 }}
194 }};
195 }
196
197 const IndicatorFunction &ind_;
198 std::size_t maxLevel_;
199 const Geometry& geometry_;
200 double volume_;
201
202 std::vector<QuadPoint> qps_;
203};
204
205} // end namespace Dumux
206
207#endif
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
Definition: adapt.hh:29
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