3.6-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
51template<class Geometry, class IndicatorFunction>
53{
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>;
58
59 class QuadPoint
60 {
61 LocalPosition localPos_;
62 double weight_;
63 std::size_t indicator_;
64
65 public:
66 QuadPoint(LocalPosition&& localPos, double weight, std::size_t i)
67 : localPos_(weight*std::move(localPos))
68 , weight_(weight)
69 , indicator_(i)
70 {}
71
72 void add(LocalPosition&& localPos, double weight)
73 {
74 localPos_ += weight*localPos;
75 weight_ += weight;
76 }
77
78 void finalize()
79 {
80 localPos_ /= weight_;
81 }
82
83 const LocalPosition& position() const { return localPos_; }
84 double weight() const { return weight_; }
85 std::size_t indicator() const { return indicator_; }
86 };
87
88 class Triangle
89 {
90 Corners corners_;
91 public:
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]; }
95
96 GlobalPosition center() const
97 {
98 GlobalPosition center(0.0);
99 for (const auto& c : corners_)
100 center += c;
101 center /= corners_.size();
102 return center;
103 }
104
105 auto volume() const
106 {
107 const auto ab = corners_[1] - corners_[0];
108 const auto ac = corners_[2] - corners_[0];
109 return crossProduct(ab, ac).two_norm()*0.5;
110 }
111 };
112
113public:
114 LocalRefinementSimplexQuadrature(const Geometry& geo, std::size_t maxLevel, const IndicatorFunction& i)
115 : ind_(i)
116 , maxLevel_(maxLevel)
117 , geometry_(geo)
118 {
119 static constexpr int dim = Geometry::mydimension;
120 static_assert(dim == 2, "Only triangles are supported so far");
121
122 if (geo.corners() != (dim+1))
123 DUNE_THROW(Dune::InvalidStateException, "Only simplex geometries are allowed");
124
125 // evaluate indicator
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();
129
130 // add new sub-qps recursively by adaptive refinement
131 addSimplex_(tri, 0, triple, triple);
132
133 for (auto& qp : qps_)
134 qp.finalize();
135 }
136
137 auto begin() const { return qps_.begin(); }
138 auto end() const { return qps_.end(); }
139 auto size() const { return qps_.size(); }
140
141private:
142 void addSimplex_(const Triangle& tri, std::size_t level, const IndicatorTriple& triple, const IndicatorTriple& childTriple)
143 {
144 // if all indicator are the same just add the triangle
145 if (std::all_of(childTriple.begin(), childTriple.end(), [a0=childTriple[0]](auto a){ return (a == a0); }))
146 {
147 // the volumes need to sum up to 0.5 (triangle reference element volume)
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_);
151 else
152 qps_.emplace_back(geometry_.local(tri.center()), 0.5*tri.volume()/volume_, childTriple[0]);
153 }
154
155 // if they are different but the maximum level is reached, split volume in three
156 else if (level == maxLevel_)
157 {
158 for (int i = 0; i < 3; ++i)
159 {
160 // the volumes need to sum up to 0.5 (triangle reference element volume)
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);
164 else
165 qps_.emplace_back(geometry_.local(tri[i]), 0.5*tri.volume()/volume_/3.0, childTriple[i]);
166 }
167 }
168
169 // otherwise refine and go recursively over all children
170 else
171 {
172 const auto children = refine(tri);
173 for (const auto& c : children)
174 {
175 // make sure to recompute indicator every time since new indices might come up
176 const auto grandChildTriple = IndicatorTriple{{ ind_(c[0]), ind_(c[1]), ind_(c[2]) }};
177 addSimplex_(c, level+1, triple, grandChildTriple);
178 }
179 }
180 }
181
182 std::array<Triangle, 4> refine(const Triangle& tri) const
183 {
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]);
187
188 return {{
189 {{ tri[0], p01, p20 }},
190 {{ tri[1], p01, p12 }},
191 {{ tri[2], p12, p20 }},
192 {{ p12, p01, p20 }}
193 }};
194 }
195
196 const IndicatorFunction &ind_;
197 std::size_t maxLevel_;
198 const Geometry& geometry_;
199 double volume_;
200
201 std::vector<QuadPoint> qps_;
202};
203
204} // end namespace Dumux
205
206#endif
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