4 #ifndef DUNE_GEOMETRY_AXISALIGNED_CUBE_GEOMETRY_HH
5 #define DUNE_GEOMETRY_AXISALIGNED_CUBE_GEOMETRY_HH
13 #include <dune/common/fvector.hh>
14 #include <dune/common/fmatrix.hh>
15 #include <dune/common/diagonalmatrix.hh>
16 #include <dune/common/unused.hh>
47 template <
class CoordType,
unsigned int dim,
unsigned int coorddim>
75 typedef typename std::conditional<dim==coorddim,
76 DiagonalMatrix<ctype,dim>,
85 typedef typename std::conditional<dim==coorddim,
86 DiagonalMatrix<ctype,dim>,
94 const Dune::FieldVector<ctype,coorddim> upper)
100 axes_ = (1<<coorddim)-1;
111 const Dune::FieldVector<ctype,coorddim> upper,
112 const std::bitset<coorddim>& axes)
117 assert(axes.count()==dim);
118 for (
size_t i=0; i<coorddim; i++)
120 upper_[i] = lower_[i];
134 lower_ = other.lower_;
135 upper_ = other.upper_;
150 if (dim == coorddim) {
151 for (
size_t i=0; i<coorddim; i++)
152 result[i] = lower_[i] +
local[i]*(upper_[i] - lower_[i]);
157 for (
size_t i=0; i<coorddim; i++)
158 result[i] = (axes_[i])
159 ? lower_[i] +
local[lc++]*(upper_[i] - lower_[i])
169 if (dim == coorddim) {
170 for (
size_t i=0; i<dim; i++)
171 result[i] = (
global[i] - lower_[i]) / (upper_[i] - lower_[i]);
172 }
else if (dim != 0) {
174 for (
size_t i=0; i<coorddim; i++)
176 result[lc++] = (
global[i] - lower_[i]) / (upper_[i] - lower_[i]);
221 for (
size_t i=0; i<coorddim; i++)
222 result[i] = 0.5 * (lower_[i] + upper_[i]);
237 if (dim == coorddim) {
238 for (
size_t i=0; i<coorddim; i++)
239 result[i] = (k & (1<<i)) ? upper_[i] : lower_[i];
243 unsigned int mask = 1;
245 for (
size_t i=0; i<coorddim; i++) {
247 result[i] = lower_[i];
249 result[i] = (k & mask) ? upper_[i] : lower_[i];
263 if (dim == coorddim) {
264 for (
size_t i=0; i<dim; i++)
265 vol *= upper_[i] - lower_[i];
267 }
else if (dim != 0) {
268 for (
size_t i=0; i<coorddim; i++)
270 vol *= upper_[i] - lower_[i];
290 for (
size_t i=0; i<dim; i++)
301 for (
size_t i=0; i<coorddim; i++)
309 for (
size_t i=0; i<dim; i++)
320 for (
size_t i=0; i<coorddim; i++)
325 Dune::FieldVector<ctype,coorddim> lower_;
327 Dune::FieldVector<ctype,coorddim> upper_;
329 std::bitset<coorddim> axes_;
A unique label for each type of element that can occur in a grid.
Definition: affinegeometry.hh:19
This class provides access to geometric and topological properties of a reference element.
Definition: referenceelements.hh:29
static const ReferenceElement< ctype, dim > & cube()
get hypercube reference elements
Definition: referenceelements.hh:472
A geometry implementation for axis-aligned hypercubes.
Definition: axisalignedcubegeometry.hh:49
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:110
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:93
AxisAlignedCubeGeometry(const Dune::FieldVector< ctype, coorddim > lower)
Constructor from a single point only.
Definition: axisalignedcubegeometry.hh:127
friend const ReferenceElement< ctype, dim > & referenceElement(const AxisAlignedCubeGeometry &geometry)
Definition: axisalignedcubegeometry.hh:281
GlobalCoordinate corner(int k) const
Return world coordinates of the k-th corner of the element.
Definition: axisalignedcubegeometry.hh:234
JacobianTransposed jacobianTransposed(DUNE_UNUSED const LocalCoordinate &local) const
Jacobian transposed of the transformation from local to global coordinates.
Definition: axisalignedcubegeometry.hh:182
@ mydimension
Definition: axisalignedcubegeometry.hh:55
FieldVector< ctype, dim > LocalCoordinate
Type used for a vector of element coordinates.
Definition: axisalignedcubegeometry.hh:64
AxisAlignedCubeGeometry & operator=(const AxisAlignedCubeGeometry &other)
Assignment operator.
Definition: axisalignedcubegeometry.hh:132
FieldVector< ctype, coorddim > GlobalCoordinate
Type used for a vector of world coordinates.
Definition: axisalignedcubegeometry.hh:67
JacobianInverseTransposed jacobianInverseTransposed(DUNE_UNUSED const LocalCoordinate &local) const
Jacobian transposed of the transformation from local to global coordinates.
Definition: axisalignedcubegeometry.hh:194
LocalCoordinate local(const GlobalCoordinate &global) const
Map a point in global (world) coordinates to element coordinates.
Definition: axisalignedcubegeometry.hh:166
CoordType ctype
Type used for single coordinate coefficients.
Definition: axisalignedcubegeometry.hh:61
@ coorddimension
Definition: axisalignedcubegeometry.hh:58
std::conditional< dim==coorddim, DiagonalMatrix< ctype, dim >, FieldMatrix< ctype, dim, coorddim > >::type JacobianTransposed
Return type of jacobianTransposed.
Definition: axisalignedcubegeometry.hh:77
GeometryType type() const
Type of the cube. Here: a hypercube of the correct dimension.
Definition: axisalignedcubegeometry.hh:141
int corners() const
Return the number of corners of the element.
Definition: axisalignedcubegeometry.hh:228
ctype integrationElement(DUNE_UNUSED const LocalCoordinate &local) const
Return the integration element, i.e., the determinant term in the integral transformation formula.
Definition: axisalignedcubegeometry.hh:208
std::conditional< dim==coorddim, DiagonalMatrix< ctype, dim >, FieldMatrix< ctype, coorddim, dim > >::type JacobianInverseTransposed
Return type of jacobianInverseTransposed.
Definition: axisalignedcubegeometry.hh:87
GlobalCoordinate center() const
Return center of mass of the element.
Definition: axisalignedcubegeometry.hh:214
ctype volume() const
Return the element volume.
Definition: axisalignedcubegeometry.hh:260
bool affine() const
Return if the element is affine. Here: yes.
Definition: axisalignedcubegeometry.hh:276
GlobalCoordinate global(const LocalCoordinate &local) const
Map a point in local (element) coordinates to world coordinates.
Definition: axisalignedcubegeometry.hh:147
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:25
@ cube
Cube element in any nonnegative dimension.
Definition: type.hh:31