version 3.11-dev
quadraturerules.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// SPDX-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_DISCRETIZATION_QUADRATURE_RULES_HH
13#define DUMUX_DISCRETIZATION_QUADRATURE_RULES_HH
14
15#include <array>
16#include <ranges>
17#include <dune/common/std/type_traits.hh>
18#include <dune/geometry/quadraturerules.hh>
19
20#include <dumux/common/tag.hh>
22#include <dumux/common/concepts/ipdata_.hh>
27
28namespace Dumux {
29
34template<class IpData>
36{
37 using Scalar = typename IpData::GlobalPosition::value_type;
38
39public:
40 QuadraturePointData(Scalar w, IpData ip)
41 : weight_(w), ipData_(std::move(ip))
42 {}
43
44 const Scalar& weight() const { return weight_; }
45 const IpData& ipData() const { return ipData_; }
46
47private:
48 Scalar weight_;
49 IpData ipData_;
50};
51
52namespace QuadratureRules {
53
58struct MidpointQuadrature : public Utility::Tag<MidpointQuadrature> {};
59
64template<int order>
65struct DuneQuadrature : public Utility::Tag<DuneQuadrature<order>>
66{
67 static constexpr int quadratureOrder = order;
68};
69
70} // end namespace QuadratureRules
71
72namespace CVFE {
73
78template<class GridView,
81 class ElementRule = QuadratureRules::MidpointQuadrature,
82 class IntersectionRule = QuadratureRules::MidpointQuadrature>
84{
85 using ScvQuadratureRule = ScvRule;
86 using ScvfQuadratureRule = ScvfRule;
87 using ElementQuadratureRule = ElementRule;
88 using IntersectionQuadratureRule = IntersectionRule;
89};
90
91namespace Detail {
92
93template<class FVElementGeometry>
94using DefinesScvQuadratureRuleType = typename FVElementGeometry::ScvQuadratureRule;
95
96template<class FVElementGeometry>
97using DefinesScvfQuadratureRuleType = typename FVElementGeometry::ScvfQuadratureRule;
98
99template<class FVElementGeometry>
100using DefinesElementQuadratureRuleType = typename FVElementGeometry::ElementQuadratureRule;
101
102template<class FVElementGeometry>
103using DefinesIntersectionQuadratureRuleType = typename FVElementGeometry::IntersectionQuadratureRule;
104
106template<class FVElementGeometry>
109 FVElementGeometry>;
110
112template<class FVElementGeometry>
115 FVElementGeometry>;
116
118template<class FVElementGeometry>
121 FVElementGeometry>;
122
124template<class FVElementGeometry>
127 FVElementGeometry>;
128
130template<class GridView>
132 typename GridView::template Codim<0>::Entity::Geometry::LocalCoordinate,
133 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate,
135 >;
136
138template<class GridView>
140 typename GridView::template Codim<0>::Entity::Geometry::LocalCoordinate,
141 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate
142 >;
143
144} // end namespace Detail
145
147template<class FVElementGeometry>
148auto quadratureRule(const FVElementGeometry& fvGeometry,
149 const typename FVElementGeometry::SubControlVolume& scv,
151{
152 using GridView = typename FVElementGeometry::GridGeometry::GridView;
154 const auto& elementGeo = fvGeometry.elementGeometry();
155 const auto ipGlobal = scv.center();
156 const auto ipLocal = elementGeo.local(ipGlobal);
157 const auto localDofIdx = scv.localDofIndex();
158
159 std::array<QuadraturePointData<QIpData>, 1> result{{
160 {scv.volume(), QIpData(0, ipLocal, ipGlobal, localDofIdx)}
161 }};
162 return result;
163}
164
166template<int order, class FVElementGeometry>
167auto quadratureRule(const FVElementGeometry& fvGeometry,
168 const typename FVElementGeometry::SubControlVolume& scv,
170{
171 using GridGeometry = typename FVElementGeometry::GridGeometry;
172 using GridView = typename GridGeometry::GridView;
175
176 auto scvGeo = fvGeometry.geometry(scv);
177 const auto& elementGeo = fvGeometry.elementGeometry();
178 auto quad = Dune::QuadratureRules<typename GridView::ctype, GridView::dimension>::rule(scvGeo.type(), order);
179 const auto localDofIdx = scv.localDofIndex();
180
181 return std::views::iota(0u, quad.size())
182 | std::views::transform([quad = std::move(quad), scvGeo = std::move(scvGeo), &elementGeo, localDofIdx](const auto idx) {
183 const auto& qp = quad[idx];
184 const auto ipGlobal = scvGeo.global(qp.position());
185 const auto ipLocal = elementGeo.local(ipGlobal);
187 qp.weight() * Extrusion::integrationElement(scvGeo, qp.position()),
188 QIpData(idx, ipLocal, ipGlobal, localDofIdx)
189 };
190 });
191}
192
194template<class FVElementGeometry>
195auto quadratureRule(const FVElementGeometry& fvGeometry,
196 const typename FVElementGeometry::SubControlVolumeFace& scvf,
198{
199 using GridView = typename FVElementGeometry::GridGeometry::GridView;
202 const auto& elementGeo = fvGeometry.elementGeometry();
203 const auto ipGlobal = scvf.ipGlobal();
204 const auto ipLocal = elementGeo.local(ipGlobal);
205
206 std::array<QuadraturePointData<QFaceIpData>, 1> result{{
207 {scvf.area(), QFaceIpData(0, scvf.unitOuterNormal(), scvf.index(), ipLocal, ipGlobal)}
208 }};
209 return result;
210}
211
213template<int order, class FVElementGeometry>
214auto quadratureRule(const FVElementGeometry& fvGeometry,
215 const typename FVElementGeometry::SubControlVolumeFace& scvf,
217{
218 using GridGeometry = typename FVElementGeometry::GridGeometry;
219 using GridView = typename GridGeometry::GridView;
222
224
225 auto scvfGeo = fvGeometry.geometry(scvf);
226 const auto& elementGeo = fvGeometry.elementGeometry();
227 auto quad = Dune::QuadratureRules<typename GridView::ctype, GridView::dimension-1>::rule(scvfGeo.type(), order);
228 auto normal = scvf.unitOuterNormal();
229 const auto scvfIndex = scvf.index();
230
231 return std::views::iota(0u, quad.size())
232 | std::views::transform([quad = std::move(quad), scvfGeo = std::move(scvfGeo), &elementGeo, normal = std::move(normal), scvfIndex](const auto idx) {
233 const auto& qp = quad[idx];
234 const auto ipGlobal = scvfGeo.global(qp.position());
235 const auto ipLocal = elementGeo.local(ipGlobal);
237 qp.weight() * Extrusion::integrationElement(scvfGeo, qp.position()),
238 QFaceIpData(idx, normal, scvfIndex, ipLocal, ipGlobal)
239 };
240 });
241}
242
244template<class FVElementGeometry>
245auto quadratureRule(const FVElementGeometry& fvGeometry,
246 const typename FVElementGeometry::Element& element,
248{
249 using GridGeometry = typename FVElementGeometry::GridGeometry;
250 using GridView = typename GridGeometry::GridView;
252
253 const auto& elementGeo = fvGeometry.elementGeometry();
254 const auto ipGlobal = elementGeo.center();
255 const auto ipLocal = elementGeo.local(ipGlobal);
256
257 std::array<Dumux::QuadraturePointData<QIpData>, 1> result{{
258 {elementGeo.volume(), QIpData(0, ipLocal, ipGlobal)}
259 }};
260 return result;
261}
262
264template<int order, class FVElementGeometry>
265auto quadratureRule(const FVElementGeometry& fvGeometry,
266 const typename FVElementGeometry::Element& element,
268{
269 using GridGeometry = typename FVElementGeometry::GridGeometry;
270 using GridView = typename GridGeometry::GridView;
273
274 const auto& elementGeo = fvGeometry.elementGeometry();
275 auto quad = Dune::QuadratureRules<typename GridView::ctype, GridView::dimension>::rule(elementGeo.type(), order);
276
277 return std::views::iota(0u, quad.size())
278 | std::views::transform([quad = std::move(quad), &elementGeo](const auto idx) {
279 const auto& qp = quad[idx];
280 const auto ipGlobal = elementGeo.global(qp.position());
281 const auto ipLocal = qp.position();
283 qp.weight() * Extrusion::integrationElement(elementGeo, qp.position()),
284 QIpData(idx, ipLocal, ipGlobal)
285 };
286 });
287}
288
290template<class FVElementGeometry>
291auto quadratureRule(const FVElementGeometry& fvGeometry,
292 const typename FVElementGeometry::GridGeometry::GridView::Intersection& is,
294{
295 // per default we don't add any intersection information to the ipData
296 using GridGeometry = typename FVElementGeometry::GridGeometry;
297 using GridView = typename GridGeometry::GridView;
300 BFlag,
303
304 const auto& elementGeo = fvGeometry.elementGeometry();
305 const auto isGeometry = is.geometry();
306 const auto ipGlobal = isGeometry.center();
307 const auto ipLocal = elementGeo.local(ipGlobal);
308
309 std::array<Dumux::QuadraturePointData<QFaceIpData>, 1> result{{
310 {isGeometry.volume(),
311 QFaceIpData(0, is.centerUnitOuterNormal(), BFlag{ is }, is.indexInInside(), ipLocal, ipGlobal)}
312 }};
313 return result;
314}
315
317template<int order, class FVElementGeometry>
318auto quadratureRule(const FVElementGeometry& fvGeometry,
319 const typename FVElementGeometry::GridGeometry::GridView::Intersection& is,
321{
322 // per default we don't add any intersection information to the ipData
323 using GridGeometry = typename FVElementGeometry::GridGeometry;
324 using GridView = typename GridGeometry::GridView;
327 BFlag,
331
332 const auto& elementGeo = fvGeometry.elementGeometry();
333 auto isGeometry = is.geometry();
334 auto quad = Dune::QuadratureRules<typename GridView::ctype, GridView::dimension-1>::rule(isGeometry.type(), order);
335 auto normal = is.centerUnitOuterNormal();
336 auto bFlag = BFlag{ is };
337 auto index = is.indexInInside();
338
339 return std::views::iota(0u, quad.size())
340 | std::views::transform(
341 [quad = std::move(quad), isGeometry = std::move(isGeometry), &elementGeo, normal = std::move(normal),
342 bFlag = std::move(bFlag), index](const auto idx) {
343 const auto& qp = quad[idx];
344 const auto ipGlobal = isGeometry.global(qp.position());
345 const auto ipLocal = elementGeo.local(ipGlobal);
347 qp.weight() * Extrusion::integrationElement(isGeometry, qp.position()),
348 QFaceIpData(idx, normal, bFlag, index, ipLocal, ipGlobal)
349 };
350 });
351}
352
354template<class FVElementGeometry>
355auto quadratureRule(const FVElementGeometry& fvGeometry, const typename FVElementGeometry::SubControlVolume& scv)
356{
357 return quadratureRule(fvGeometry, scv, typename Detail::ScvQuadratureRuleType<FVElementGeometry>{});
358}
359
361template<class FVElementGeometry>
362auto quadratureRule(const FVElementGeometry& fvGeometry, const typename FVElementGeometry::SubControlVolumeFace& scvf)
363{
364 return quadratureRule(fvGeometry, scvf, typename Detail::ScvfQuadratureRuleType<FVElementGeometry>{});
365}
366
368template<class FVElementGeometry>
369auto quadratureRule(const FVElementGeometry& fvGeometry, const typename FVElementGeometry::Element& element)
370{
371 return quadratureRule(fvGeometry, element, typename Detail::ElementQuadratureRuleType<FVElementGeometry>{});
372}
373
375template<class FVElementGeometry>
376auto quadratureRule(const FVElementGeometry& fvGeometry, const typename FVElementGeometry::GridGeometry::GridView::Intersection& intersection)
377{
378 return quadratureRule(fvGeometry, intersection, typename Detail::IntersectionQuadratureRuleType<FVElementGeometry>{});
379}
380
381} // end namespace CVFE
382
383} // end namespace Dumux
384
385#endif
Boundary flag to store e.g. in sub control volume faces.
Boundary flag to store e.g. in sub control volume faces.
Definition: boundaryflag.hh:58
An interpolation point related to a face of an element.
Definition: cvfe/interpolationpointdata.hh:103
Wraps interpolation point data and adds a quadrature point index for use in quadrature loops.
Definition: cvfe/interpolationpointdata.hh:135
An interpolation point related to an element that includes global and local positions.
Definition: cvfe/interpolationpointdata.hh:25
An interpolation point related to a localDof of an element, giving its global and local positions.
Definition: cvfe/interpolationpointdata.hh:54
Traits extracting the public Extrusion type from T Defaults to NoExtrusion if no such type is found.
Definition: extrusion.hh:155
An interpolation point related to an intersection.
Definition: fem/interpolationpointdata.hh:107
Quadrature point data with weight and interpolation point info.
Definition: quadraturerules.hh:36
const Scalar & weight() const
Definition: quadraturerules.hh:44
QuadraturePointData(Scalar w, IpData ip)
Definition: quadraturerules.hh:40
const IpData & ipData() const
Definition: quadraturerules.hh:45
Classes representing interpolation point data for control-volume finite element schemes.
Helper classes to compute the integration elements.
Shape functions and gradients at an interpolation point.
Vector normal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: normal.hh:26
Defines the index types used for grid and local indices.
typename FVElementGeometry::ScvQuadratureRule DefinesScvQuadratureRuleType
Definition: quadraturerules.hh:94
Dune::Std::detected_or_t< QuadratureRules::MidpointQuadrature, DefinesScvfQuadratureRuleType, FVElementGeometry > ScvfQuadratureRuleType
Get ScvfQuadratureRule type (default is MidpointQuadrature)
Definition: quadraturerules.hh:115
typename FVElementGeometry::ElementQuadratureRule DefinesElementQuadratureRuleType
Definition: quadraturerules.hh:100
Dune::Std::detected_or_t< QuadratureRules::MidpointQuadrature, DefinesElementQuadratureRuleType, FVElementGeometry > ElementQuadratureRuleType
Get ElementQuadratureRule type (default is MidpointQuadrature)
Definition: quadraturerules.hh:121
Dune::Std::detected_or_t< QuadratureRules::MidpointQuadrature, DefinesScvQuadratureRuleType, FVElementGeometry > ScvQuadratureRuleType
Get ScvQuadratureRule type (default is MidpointQuadrature)
Definition: quadraturerules.hh:109
Dune::Std::detected_or_t< QuadratureRules::MidpointQuadrature, DefinesIntersectionQuadratureRuleType, FVElementGeometry > IntersectionQuadratureRuleType
Get IntersectionQuadratureRule type (default is MidpointQuadrature)
Definition: quadraturerules.hh:127
typename FVElementGeometry::ScvfQuadratureRule DefinesScvfQuadratureRuleType
Definition: quadraturerules.hh:97
typename FVElementGeometry::IntersectionQuadratureRule DefinesIntersectionQuadratureRuleType
Definition: quadraturerules.hh:103
auto quadratureRule(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolume &scv, QuadratureRules::MidpointQuadrature)
Midpoint quadrature for scv.
Definition: quadraturerules.hh:148
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
Quadrature rule traits for discretization schemes.
Definition: quadraturerules.hh:84
IntersectionRule IntersectionQuadratureRule
Definition: quadraturerules.hh:88
ScvfRule ScvfQuadratureRule
Definition: quadraturerules.hh:86
ScvRule ScvQuadratureRule
Definition: quadraturerules.hh:85
ElementRule ElementQuadratureRule
Definition: quadraturerules.hh:87
unsigned int LocalIndex
Definition: indextraits.hh:28
Dune-based quadrature rule with specified order.
Definition: quadraturerules.hh:66
static constexpr int quadratureOrder
Definition: quadraturerules.hh:67
Midpoint quadrature rule that uses scv/scvf centers.
Definition: quadraturerules.hh:58
Helper class to create (named and comparable) tagged types Tags any given type. The tagged type is eq...
Definition: tag.hh:30
Helper class to create (named and comparable) tagged types.