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>
48 template <
class CoordType,
unsigned int dim,
unsigned int coorddim>
76 typedef typename conditional<dim==coorddim,
77 DiagonalMatrix<ctype,dim>,
86 typedef typename conditional<dim==coorddim,
87 DiagonalMatrix<ctype,dim>,
95 const Dune::FieldVector<ctype,coorddim> upper)
101 axes_ = (1<<coorddim)-1;
112 const Dune::FieldVector<ctype,coorddim> upper,
113 const std::bitset<coorddim>& axes)
118 assert(axes.count()==dim);
119 for (
size_t i=0; i<coorddim; i++)
121 upper_[i] = lower_[i];
135 lower_ = other.lower_;
136 upper_ = other.upper_;
151 if (dim == coorddim) {
152 for (
size_t i=0; i<coorddim; i++)
153 result[i] = lower_[i] + local[i]*(upper_[i] - lower_[i]);
158 for (
size_t i=0; i<coorddim; i++)
159 result[i] = (axes_[i])
160 ? lower_[i] + local[lc++]*(upper_[i] - lower_[i])
170 if (dim == coorddim) {
171 for (
size_t i=0; i<dim; i++)
172 result[i] = (global[i] - lower_[i]) / (upper_[i] - lower_[i]);
173 }
else if (dim != 0) {
175 for (
size_t i=0; i<coorddim; i++)
177 result[lc++] = (global[i] - lower_[i]) / (upper_[i] - lower_[i]);
222 for (
size_t i=0; i<coorddim; i++)
223 result[i] = 0.5 * (lower_[i] + upper_[i]);
238 if (dim == coorddim) {
239 for (
size_t i=0; i<coorddim; i++)
240 result[i] = (k & (1<<i)) ? upper_[i] : lower_[i];
244 unsigned int mask = 1;
246 for (
size_t i=0; i<coorddim; i++) {
248 result[i] = lower_[i];
250 result[i] = (k & mask) ? upper_[i] : lower_[i];
264 if (dim == coorddim) {
265 for (
size_t i=0; i<dim; i++)
266 vol *= upper_[i] - lower_[i];
268 }
else if (dim != 0) {
269 for (
size_t i=0; i<coorddim; i++)
271 vol *= upper_[i] - lower_[i];
286 for (
size_t i=0; i<dim; i++)
287 jacobianTransposed.diagonal()[i] = upper_[i] - lower_[i];
297 for (
size_t i=0; i<coorddim; i++)
299 jacobianTransposed[lc++][i] = upper_[i] - lower_[i];
305 for (
size_t i=0; i<dim; i++)
306 jacobianInverseTransposed.diagonal()[i] = 1.0 / (upper_[i] - lower_[i]);
316 for (
size_t i=0; i<coorddim; i++)
318 jacobianInverseTransposed[i][lc++] = 1.0 / (upper_[i] - lower_[i]);
321 Dune::FieldVector<ctype,coorddim> lower_;
323 Dune::FieldVector<ctype,coorddim> upper_;
325 std::bitset<coorddim> axes_;
JacobianInverseTransposed jacobianInverseTransposed(DUNE_UNUSED const LocalCoordinate &local) const
Jacobian transposed of the transformation from local to global coordinates.
Definition: axisalignedcubegeometry.hh:195
conditional< dim==coorddim, DiagonalMatrix< ctype, dim >, FieldMatrix< ctype, coorddim, dim > >::type JacobianInverseTransposed
Return type of jacobianInverseTransposed.
Definition: axisalignedcubegeometry.hh:88
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:111
ctype volume() const
Return the element volume.
Definition: axisalignedcubegeometry.hh:261
FieldVector< ctype, coorddim > GlobalCoordinate
Type used for a vector of world coordinates.
Definition: axisalignedcubegeometry.hh:68
AxisAlignedCubeGeometry(const Dune::FieldVector< ctype, coorddim > lower)
Constructor from a single point only.
Definition: axisalignedcubegeometry.hh:128
Definition: axisalignedcubegeometry.hh:59
GlobalCoordinate corner(int k) const
Return world coordinates of the k-th corner of the element.
Definition: axisalignedcubegeometry.hh:235
AxisAlignedCubeGeometry & operator=(const AxisAlignedCubeGeometry &other)
Assignment operator.
Definition: axisalignedcubegeometry.hh:133
conditional< dim==coorddim, DiagonalMatrix< ctype, dim >, FieldMatrix< ctype, dim, coorddim > >::type JacobianTransposed
Return type of jacobianTransposed.
Definition: axisalignedcubegeometry.hh:78
GlobalCoordinate center() const
Return center of mass of the element.
Definition: axisalignedcubegeometry.hh:215
Cube element in any nonnegative dimension.
Definition: type.hh:31
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:94
FieldVector< ctype, dim > LocalCoordinate
Type used for a vector of element coordinates.
Definition: axisalignedcubegeometry.hh:65
GlobalCoordinate global(const LocalCoordinate &local) const
Map a point in local (element) coordinates to world coordinates.
Definition: axisalignedcubegeometry.hh:148
A unique label for each type of element that can occur in a grid.
JacobianTransposed jacobianTransposed(DUNE_UNUSED const LocalCoordinate &local) const
Jacobian transposed of the transformation from local to global coordinates.
Definition: axisalignedcubegeometry.hh:183
int corners() const
Return the number of corners of the element.
Definition: axisalignedcubegeometry.hh:229
GeometryType type() const
Type of the cube. Here: a hypercube of the correct dimension.
Definition: axisalignedcubegeometry.hh:142
LocalCoordinate local(const GlobalCoordinate &global) const
Map a point in global (world) coordinates to element coordinates.
Definition: axisalignedcubegeometry.hh:167
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:24
bool affine() const
Return if the element is affine. Here: yes.
Definition: axisalignedcubegeometry.hh:277
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:209
Definition: axisalignedcubegeometry.hh:56
CoordType ctype
Type used for single coordinate coefficients.
Definition: axisalignedcubegeometry.hh:62
A geometry implementation for axis-aligned hypercubes.
Definition: axisalignedcubegeometry.hh:49