12#ifndef DUMUX_GEOMETRY_VOLUME_HH
13#define DUMUX_GEOMETRY_VOLUME_HH
19#include <dune/common/exceptions.hh>
20#include <dune/geometry/type.hh>
21#include <dune/geometry/quadraturerules.hh>
40template<
int dim,
class CornerF>
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");
48 if constexpr (coordDim > 3)
49 return std::numeric_limits<ctype>::quiet_NaN();
51 if constexpr (dim == 0)
54 else if constexpr (dim == 1)
55 return (c(1)-c(0)).two_norm();
57 else if constexpr (dim == 2)
59 if (type == Dune::GeometryTypes::triangle)
61 if constexpr (coordDim == 2)
73 else if (type == Dune::GeometryTypes::quadrilateral)
75 if constexpr (coordDim == 2)
88 return std::numeric_limits<ctype>::quiet_NaN();
91 else if constexpr (dim == 3)
93 if (type == Dune::GeometryTypes::tetrahedron)
100 else if (type == Dune::GeometryTypes::hexahedron)
103 const auto v = c(7)-c(0);
111 else if (type == Dune::GeometryTypes::pyramid)
118 return 1.0/6.0 * abs(
122 else if (type == Dune::GeometryTypes::prism)
132 return std::numeric_limits<ctype>::quiet_NaN();
135 return std::numeric_limits<ctype>::quiet_NaN();
142template<
class Geometry>
145 const auto v = convexPolytopeVolume<Geometry::mydimension>(
146 geo.type(), [&](
unsigned int i){ return geo.corner(i); }
151 return std::isnan(v) ? geo.volume() : v;
158template<
class Geometry>
159auto volume(
const Geometry& geo,
unsigned int integrationOrder = 4)
161 using ctype =
typename Geometry::ctype;
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();
174template<
class Geometry,
class Transformation>
175auto volume(
const Geometry& geo, Transformation transformation,
unsigned int integrationOrder = 4)
177 using ctype =
typename Geometry::ctype;
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();
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:671
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:700
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.