version 3.8
volume.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-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_GEOMETRY_VOLUME_HH
13#define DUMUX_GEOMETRY_VOLUME_HH
14
15#include <cmath>
16#include <limits>
17#include <type_traits>
18
19#include <dune/common/exceptions.hh>
20#include <dune/geometry/type.hh>
21#include <dune/geometry/quadraturerules.hh>
22
23#include <dumux/common/math.hh>
24
25namespace Dumux {
26
40template<int dim, class CornerF>
41auto convexPolytopeVolume(Dune::GeometryType type, const CornerF& c)
42{
43 using ctype = typename std::decay_t<decltype(c(0))>::value_type;
44 static constexpr int coordDim = std::decay_t<decltype(c(0))>::dimension;
45 static_assert(coordDim >= dim, "Coordinate dimension has to be larger than geometry dimension");
46
47 // not implemented for coordinate dimension larger than 3
48 if constexpr (coordDim > 3)
49 return std::numeric_limits<ctype>::quiet_NaN();
50
51 if constexpr (dim == 0)
52 return 1.0;
53
54 else if constexpr (dim == 1)
55 return (c(1)-c(0)).two_norm();
56
57 else if constexpr (dim == 2)
58 {
59 if (type == Dune::GeometryTypes::triangle)
60 {
61 if constexpr (coordDim == 2)
62 {
63 // make sure we are using positive volumes
64 // the cross product of edge vectors might be negative,
65 // depending on the element orientation
66 using std::abs;
67 return 0.5*abs(Dumux::crossProduct(c(1)-c(0), c(2)-c(0)));
68 }
69 else // coordDim == 3
70 return 0.5*Dumux::crossProduct(c(1)-c(0), c(2)-c(0)).two_norm();
71
72 }
73 else if (type == Dune::GeometryTypes::quadrilateral)
74 {
75 if constexpr (coordDim == 2)
76 {
77 // make sure we are using positive volumes
78 // the cross product of diagonals might be negative,
79 // depending on the element orientation
80 using std::abs;
81 return 0.5*abs(Dumux::crossProduct(c(3)-c(0), c(2)-c(1)));
82 }
83 else // coordDim == 3
84 return 0.5*Dumux::crossProduct(c(3)-c(0), c(2)-c(1)).two_norm();
85
86 }
87 else
88 return std::numeric_limits<ctype>::quiet_NaN();
89 }
90
91 else if constexpr (dim == 3)
92 {
93 if (type == Dune::GeometryTypes::tetrahedron)
94 {
95 using std::abs;
96 return 1.0/6.0 * abs(
97 Dumux::tripleProduct(c(3)-c(0), c(1)-c(0), c(2)-c(0))
98 );
99 }
100 else if (type == Dune::GeometryTypes::hexahedron)
101 {
102 // after Grandy 1997, Efficient computation of volume of hexahedron
103 const auto v = c(7)-c(0);
104 using std::abs;
105 return 1.0/6.0 * (
106 abs(Dumux::tripleProduct(v, c(1)-c(0), c(3)-c(5)))
107 + abs(Dumux::tripleProduct(v, c(4)-c(0), c(5)-c(6)))
108 + abs(Dumux::tripleProduct(v, c(2)-c(0), c(6)-c(3)))
109 );
110 }
111 else if (type == Dune::GeometryTypes::pyramid)
112 {
113 // 1/3 * base * height
114 // for base see case Dune::GeometryTypes::quadrilateral above
115 // = 1/3 * (1/2 * norm(ADxBC)) * ((ADxBC)/norm(AD x BC) ⋅ AE)
116 // = 1/6 * (AD x BC) ⋅ AE
117 using std::abs;
118 return 1.0/6.0 * abs(
119 Dumux::tripleProduct(c(3)-c(0), c(2)-c(1), c(4)-c(0))
120 );
121 }
122 else if (type == Dune::GeometryTypes::prism)
123 {
124 // compute as sum of a pyramid (0-1-3-4-5) and a tetrahedron (2-0-1-5)
125 using std::abs;
126 return 1.0/6.0 * (
127 abs(Dumux::tripleProduct(c(3)-c(1), c(4)-c(0), c(5)-c(0)))
128 + abs(Dumux::tripleProduct(c(5)-c(2), c(0)-c(2), c(1)-c(2)))
129 );
130 }
131 else
132 return std::numeric_limits<ctype>::quiet_NaN();
133 }
134 else
135 return std::numeric_limits<ctype>::quiet_NaN();
136}
137
142template<class Geometry>
143auto convexPolytopeVolume(const Geometry& geo)
144{
145 const auto v = convexPolytopeVolume<Geometry::mydimension>(
146 geo.type(), [&](unsigned int i){ return geo.corner(i); }
147 );
148
149 // fall back to the method of the geometry if no specialized
150 // volume function is implemented for the geometry type
151 return std::isnan(v) ? geo.volume() : v;
152}
153
158template<class Geometry>
159auto volume(const Geometry& geo, unsigned int integrationOrder = 4)
160{
161 using ctype = typename Geometry::ctype;
162 ctype volume = 0.0;
163 const auto rule = Dune::QuadratureRules<ctype, Geometry::mydimension>::rule(geo.type(), integrationOrder);
164 for (const auto& qp : rule)
165 volume += geo.integrationElement(qp.position())*qp.weight();
166 return volume;
167}
168
174template<class Geometry, class Transformation>
175auto volume(const Geometry& geo, Transformation transformation, unsigned int integrationOrder = 4)
176{
177 using ctype = typename Geometry::ctype;
178 ctype volume = 0.0;
179 const auto rule = Dune::QuadratureRules<ctype, Geometry::mydimension>::rule(geo.type(), integrationOrder);
180 for (const auto& qp : rule)
181 volume += transformation.integrationElement(geo, qp.position())*qp.weight();
182 return volume;
183}
184
185} // end namespace Dumux
186
187#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:642
Scalar tripleProduct(const Dune::FieldVector< Scalar, 3 > &vec1, const Dune::FieldVector< Scalar, 3 > &vec2, const Dune::FieldVector< Scalar, 3 > &vec3)
Triple product of three vectors in three-dimensional Euclidean space retuning scalar.
Definition: math.hh:671
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:159
auto convexPolytopeVolume(Dune::GeometryType type, const CornerF &c)
Compute the volume of several common geometry types.
Definition: volume.hh:41
Define some often used mathematical functions.
Definition: adapt.hh:17