6#ifndef DUNE_GEOMETRY_AXISALIGNED_CUBE_GEOMETRY_HH
7#define DUNE_GEOMETRY_AXISALIGNED_CUBE_GEOMETRY_HH
15#include <dune/common/fvector.hh>
16#include <dune/common/fmatrix.hh>
17#include <dune/common/diagonalmatrix.hh>
48 template <
class CoordType,
unsigned int dim,
unsigned int coorddim>
79 typedef typename std::conditional<dim==coorddim,
80 DiagonalMatrix<ctype,dim>,
89 typedef typename std::conditional<dim==coorddim,
90 DiagonalMatrix<ctype,dim>,
100 using Jacobian = std::conditional_t<dim==coorddim, DiagonalMatrix<ctype,dim>, FieldMatrix<ctype,coorddim,dim> >;
109 using JacobianInverse = std::conditional_t<dim==coorddim, DiagonalMatrix<ctype,dim>, FieldMatrix<ctype,dim,coorddim> >;
125 const Dune::FieldVector<ctype,coorddim> upper)
130 static_assert(dim==coorddim,
"Use this constructor only if dim==coorddim!");
132 axes_ = (1<<coorddim)-1;
143 const Dune::FieldVector<ctype,coorddim> upper,
144 const std::bitset<coorddim>& axes)
149 assert(axes.count()==dim);
150 for (
size_t i=0; i<coorddim; i++)
152 upper_[i] = lower_[i];
166 return GeometryTypes::cube(dim);
173 if (dim == coorddim) {
174 for (
size_t i=0; i<coorddim; i++)
175 result[i] = lower_[i] +
local[i]*(upper_[i] - lower_[i]);
176 }
else if (dim == 0) {
180 for (
size_t i=0; i<coorddim; i++)
181 result[i] = (axes_[i])
182 ? lower_[i] +
local[lc++]*(upper_[i] - lower_[i])
192 if (dim == coorddim) {
193 for (
size_t i=0; i<dim; i++)
194 result[i] = (
global[i] - lower_[i]) / (upper_[i] - lower_[i]);
195 }
else if (dim != 0) {
197 for (
size_t i=0; i<coorddim; i++)
199 result[lc++] = (
global[i] - lower_[i]) / (upper_[i] - lower_[i]);
256 for (
size_t i=0; i<coorddim; i++)
257 result[i] = CoordType(0.5) * (lower_[i] + upper_[i]);
272 if (dim == coorddim) {
273 for (
size_t i=0; i<coorddim; i++)
274 result[i] = (k & (1<<i)) ? upper_[i] : lower_[i];
275 }
else if (dim == 0) {
278 unsigned int mask = 1;
280 for (
size_t i=0; i<coorddim; i++) {
282 result[i] = lower_[i];
284 result[i] = (k & mask) ? upper_[i] : lower_[i];
298 if (dim == coorddim) {
299 for (
size_t i=0; i<dim; i++)
300 vol *= upper_[i] - lower_[i];
302 }
else if (dim != 0) {
303 for (
size_t i=0; i<coorddim; i++)
305 vol *= upper_[i] - lower_[i];
325 for (
size_t i=0; i<dim; i++)
336 for (
size_t i=0; i<coorddim; i++)
344 for (
size_t i=0; i<dim; i++)
355 for (
size_t i=0; i<coorddim; i++)
360 Dune::FieldVector<ctype,coorddim> lower_;
362 Dune::FieldVector<ctype,coorddim> upper_;
364 std::bitset<coorddim> axes_;
A unique label for each type of element that can occur in a grid.
Definition affinegeometry.hh:21
typename Container::ReferenceElement ReferenceElement
The reference element type.
Definition referenceelements.hh:146
static const ReferenceElement & cube()
get hypercube reference elements
Definition referenceelements.hh:168
A geometry implementation for axis-aligned hypercubes.
Definition axisalignedcubegeometry.hh:50
Volume volume() const
Return the element volume.
Definition axisalignedcubegeometry.hh:295
AxisAlignedCubeGeometry(const Dune::FieldVector< ctype, coorddim > lower, const Dune::FieldVector< ctype, coorddim > upper, const std::bitset< coorddim > &axes)
Constructor from a lower left and an upper right corner.
Definition axisalignedcubegeometry.hh:142
AxisAlignedCubeGeometry(const Dune::FieldVector< ctype, coorddim > lower, const Dune::FieldVector< ctype, coorddim > upper)
Constructor from a lower left and an upper right corner.
Definition axisalignedcubegeometry.hh:124
JacobianInverse jacobianInverse(const LocalCoordinate &local) const
Inverse Jacobian of the transformation from local to global coordinates.
Definition axisalignedcubegeometry.hh:235
JacobianInverseTransposed jacobianInverseTransposed(const LocalCoordinate &local) const
Inverse Jacobian transposed of the transformation from local to global coordinates.
Definition axisalignedcubegeometry.hh:217
std::conditional< dim==coorddim, DiagonalMatrix< ctype, dim >, FieldMatrix< ctype, dim, coorddim > >::type JacobianTransposed
Return type of jacobianTransposed.
Definition axisalignedcubegeometry.hh:81
AxisAlignedCubeGeometry(const Dune::FieldVector< ctype, coorddim > lower)
Constructor from a single point only.
Definition axisalignedcubegeometry.hh:159
static constexpr int mydimension
Dimension of the cube element.
Definition axisalignedcubegeometry.hh:56
friend Dune::ReferenceElements< ctype, dim >::ReferenceElement referenceElement(const AxisAlignedCubeGeometry &)
Definition axisalignedcubegeometry.hh:316
static constexpr int coorddimension
Dimension of the world space that the cube element is embedded in.
Definition axisalignedcubegeometry.hh:59
GlobalCoordinate corner(int k) const
Return world coordinates of the k-th corner of the element.
Definition axisalignedcubegeometry.hh:269
ctype Volume
Type used for volume.
Definition axisalignedcubegeometry.hh:71
FieldVector< ctype, dim > LocalCoordinate
Type used for a vector of element coordinates.
Definition axisalignedcubegeometry.hh:65
JacobianTransposed jacobianTransposed(const LocalCoordinate &local) const
Jacobian transposed of the transformation from local to global coordinates.
Definition axisalignedcubegeometry.hh:205
FieldVector< ctype, coorddim > GlobalCoordinate
Type used for a vector of world coordinates.
Definition axisalignedcubegeometry.hh:68
LocalCoordinate local(const GlobalCoordinate &global) const
Map a point in global (world) coordinates to element coordinates.
Definition axisalignedcubegeometry.hh:189
CoordType ctype
Type used for single coordinate coefficients.
Definition axisalignedcubegeometry.hh:62
std::conditional_t< dim==coorddim, DiagonalMatrix< ctype, dim >, FieldMatrix< ctype, coorddim, dim > > Jacobian
Return type of jacobian.
Definition axisalignedcubegeometry.hh:100
GeometryType type() const
Type of the cube. Here: a hypercube of the correct dimension.
Definition axisalignedcubegeometry.hh:164
int corners() const
Return the number of corners of the element.
Definition axisalignedcubegeometry.hh:263
Jacobian jacobian(const LocalCoordinate &local) const
Jacobian of the transformation from local to global coordinates.
Definition axisalignedcubegeometry.hh:229
std::conditional< dim==coorddim, DiagonalMatrix< ctype, dim >, FieldMatrix< ctype, coorddim, dim > >::type JacobianInverseTransposed
Return type of jacobianInverseTransposed.
Definition axisalignedcubegeometry.hh:91
Volume integrationElement(const LocalCoordinate &local) const
Return the integration element, i.e., the determinant term in the integral transformation formula.
Definition axisalignedcubegeometry.hh:243
GlobalCoordinate center() const
Return center of mass of the element.
Definition axisalignedcubegeometry.hh:249
AxisAlignedCubeGeometry()=default
Constructs an empty geometry.
bool affine() const
Return if the element is affine. Here: yes.
Definition axisalignedcubegeometry.hh:311
std::conditional_t< dim==coorddim, DiagonalMatrix< ctype, dim >, FieldMatrix< ctype, dim, coorddim > > JacobianInverse
Return type of jacobianInverse.
Definition axisalignedcubegeometry.hh:109
GlobalCoordinate global(const LocalCoordinate &local) const
Map a point in local (element) coordinates to world coordinates.
Definition axisalignedcubegeometry.hh:170
Unique label for each type of entities that can occur in DUNE grids.
Definition type.hh:114