3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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 *****************************************************************************/
24#ifndef DUMUX_GEOMETRY_VOLUME_HH
25#define DUMUX_GEOMETRY_VOLUME_HH
26
27#include <cmath>
28#include <limits>
29#include <type_traits>
30
31#include <dune/common/exceptions.hh>
32#include <dune/geometry/type.hh>
33#include <dune/geometry/quadraturerules.hh>
34
35#include <dumux/common/math.hh>
36
37namespace Dumux {
38
52template<int dim, class CornerF>
53auto convexPolytopeVolume(Dune::GeometryType type, const CornerF& c)
54{
55 using ctype = typename std::decay_t<decltype(c(0))>::value_type;
56 static constexpr int coordDim = std::decay_t<decltype(c(0))>::dimension;
57 static_assert(coordDim >= dim, "Coordinate dimension has to be larger than geometry dimension");
58
59 // not implemented for coordinate dimension larger than 3
60 if constexpr (coordDim > 3)
61 return std::numeric_limits<ctype>::quiet_NaN();
62
63 if constexpr (dim == 0)
64 return 1.0;
65
66 else if constexpr (dim == 1)
67 return (c(1)-c(0)).two_norm();
68
69 else if constexpr (dim == 2)
70 {
71 if (type == Dune::GeometryTypes::triangle)
72 {
73 if constexpr (coordDim == 2)
74 {
75 // make sure we are using positive volumes
76 // the cross product of edge vectors might be negative,
77 // depending on the element orientation
78 using std::abs;
79 return 0.5*abs(Dumux::crossProduct(c(1)-c(0), c(2)-c(0)));
80 }
81 else // coordDim == 3
82 return 0.5*Dumux::crossProduct(c(1)-c(0), c(2)-c(0)).two_norm();
83
84 }
85 else if (type == Dune::GeometryTypes::quadrilateral)
86 {
87 if constexpr (coordDim == 2)
88 {
89 // make sure we are using positive volumes
90 // the cross product of diagonals might be negative,
91 // depending on the element orientation
92 using std::abs;
93 return 0.5*abs(Dumux::crossProduct(c(3)-c(0), c(2)-c(1)));
94 }
95 else // coordDim == 3
96 return 0.5*Dumux::crossProduct(c(3)-c(0), c(2)-c(1)).two_norm();
97
98 }
99 else
100 return std::numeric_limits<ctype>::quiet_NaN();
101 }
102
103 else if constexpr (dim == 3)
104 {
105 if (type == Dune::GeometryTypes::tetrahedron)
106 {
107 using std::abs;
108 return 1.0/6.0 * abs(
109 Dumux::tripleProduct(c(3)-c(0), c(1)-c(0), c(2)-c(0))
110 );
111 }
112 else if (type == Dune::GeometryTypes::hexahedron)
113 {
114 // after Grandy 1997, Efficient computation of volume of hexahedron
115 const auto v = c(7)-c(0);
116 using std::abs;
117 return 1.0/6.0 * (
118 abs(Dumux::tripleProduct(v, c(1)-c(0), c(3)-c(5)))
119 + abs(Dumux::tripleProduct(v, c(4)-c(0), c(5)-c(6)))
120 + abs(Dumux::tripleProduct(v, c(2)-c(0), c(6)-c(3)))
121 );
122 }
123 else if (type == Dune::GeometryTypes::pyramid)
124 {
125 // 1/3 * base * height
126 // for base see case Dune::GeometryTypes::quadrilateral above
127 // = 1/3 * (1/2 * norm(ADxBC)) * ((ADxBC)/norm(AD x BC) ⋅ AE)
128 // = 1/6 * (AD x BC) ⋅ AE
129 using std::abs;
130 return 1.0/6.0 * abs(
131 Dumux::tripleProduct(c(3)-c(0), c(2)-c(1), c(4)-c(0))
132 );
133 }
134 else if (type == Dune::GeometryTypes::prism)
135 {
136 // compute as sum of a pyramid (0-1-3-4-5) and a tetrahedron (2-0-1-5)
137 using std::abs;
138 return 1.0/6.0 * (
139 abs(Dumux::tripleProduct(c(3)-c(1), c(4)-c(0), c(5)-c(0)))
140 + abs(Dumux::tripleProduct(c(5)-c(2), c(0)-c(2), c(1)-c(2)))
141 );
142 }
143 else
144 return std::numeric_limits<ctype>::quiet_NaN();
145 }
146 else
147 return std::numeric_limits<ctype>::quiet_NaN();
148}
149
154template<class Geometry>
155auto convexPolytopeVolume(const Geometry& geo)
156{
157 const auto v = convexPolytopeVolume<Geometry::mydimension>(
158 geo.type(), [&](unsigned int i){ return geo.corner(i); }
159 );
160
161 // fall back to the method of the geometry if no specialized
162 // volume function is implemented for the geometry type
163 return std::isnan(v) ? geo.volume() : v;
164}
165
170template<class Geometry>
171auto volume(const Geometry& geo, unsigned int integrationOrder = 4)
172{
173 using ctype = typename Geometry::ctype;
174 ctype volume = 0.0;
175 const auto rule = Dune::QuadratureRules<ctype, Geometry::mydimension>::rule(geo.type(), integrationOrder);
176 for (const auto& qp : rule)
177 volume += geo.integrationElement(qp.position())*qp.weight();
178 return volume;
179}
180
186template<class Geometry, class Transformation>
187auto volume(const Geometry& geo, Transformation transformation, unsigned int integrationOrder = 4)
188{
189 using ctype = typename Geometry::ctype;
190 ctype volume = 0.0;
191 const auto rule = Dune::QuadratureRules<ctype, Geometry::mydimension>::rule(geo.type(), integrationOrder);
192 for (const auto& qp : rule)
193 volume += transformation.integrationElement(geo, qp.position())*qp.weight();
194 return volume;
195}
196
197} // end namespace Dumux
198
199#endif
Define some often used mathematical functions.
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:171
auto convexPolytopeVolume(Dune::GeometryType type, const CornerF &c)
Compute the volume of several common geometry types.
Definition: volume.hh:53
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
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:683
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29