24#ifndef DUMUX_GEOMETRY_VOLUME_HH
25#define DUMUX_GEOMETRY_VOLUME_HH
31#include <dune/common/exceptions.hh>
32#include <dune/geometry/type.hh>
33#include <dune/geometry/quadraturerules.hh>
52template<
int dim,
class CornerF>
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");
60 if constexpr (coordDim > 3)
61 return std::numeric_limits<ctype>::quiet_NaN();
63 if constexpr (dim == 0)
66 else if constexpr (dim == 1)
67 return (c(1)-c(0)).two_norm();
69 else if constexpr (dim == 2)
71 if (type == Dune::GeometryTypes::triangle)
73 if constexpr (coordDim == 2)
85 else if (type == Dune::GeometryTypes::quadrilateral)
87 if constexpr (coordDim == 2)
100 return std::numeric_limits<ctype>::quiet_NaN();
103 else if constexpr (dim == 3)
105 if (type == Dune::GeometryTypes::tetrahedron)
108 return 1.0/6.0 * abs(
112 else if (type == Dune::GeometryTypes::hexahedron)
115 const auto v = c(7)-c(0);
123 else if (type == Dune::GeometryTypes::pyramid)
130 return 1.0/6.0 * abs(
134 else if (type == Dune::GeometryTypes::prism)
144 return std::numeric_limits<ctype>::quiet_NaN();
147 return std::numeric_limits<ctype>::quiet_NaN();
154template<
class Geometry>
157 const auto v = convexPolytopeVolume<Geometry::mydimension>(
158 geo.type(), [&](
unsigned int i){ return geo.corner(i); }
163 return std::isnan(v) ? geo.volume() : v;
170template<
class Geometry>
171auto volume(
const Geometry& geo,
unsigned int integrationOrder = 4)
173 using ctype =
typename Geometry::ctype;
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();
186template<
class Geometry,
class Transformation>
187auto volume(
const Geometry& geo, Transformation transformation,
unsigned int integrationOrder = 4)
189 using ctype =
typename Geometry::ctype;
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();
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