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)
98 axes_((1<<coorddim)-1),
99 jacobianTransposed_(0),
100 jacobianInverseTransposed_(0)
111 const Dune::FieldVector<ctype,coorddim> upper,
112 const std::bitset<coorddim>& axes)
116 jacobianTransposed_(0),
117 jacobianInverseTransposed_(0)
119 assert(axes.count()==dim);
120 for (
size_t i=0; i<coorddim; i++)
122 upper_[i] = lower_[i];
131 jacobianTransposed_(0),
132 jacobianInverseTransposed_(0)
138 lower_ = other.lower_;
139 upper_ = other.upper_;
141 jacobianTransposed_ = other.jacobianTransposed_;
142 jacobianInverseTransposed_ = other.jacobianInverseTransposed_;
156 if (dim == coorddim) {
157 for (
size_t i=0; i<coorddim; i++)
158 result[i] = lower_[i] + local[i]*(upper_[i] - lower_[i]);
163 for (
size_t i=0; i<coorddim; i++)
164 result[i] = (axes_[i])
165 ? lower_[i] + local[lc++]*(upper_[i] - lower_[i])
175 if (dim == coorddim) {
176 for (
size_t i=0; i<dim; i++)
177 result[i] = (global[i] - lower_[i]) / (upper_[i] - lower_[i]);
178 }
else if (dim != 0) {
180 for (
size_t i=0; i<coorddim; i++)
182 result[lc++] = (global[i] - lower_[i]) / (upper_[i] - lower_[i]);
191 return jacobianTransposed_;
198 return jacobianInverseTransposed_;
217 for (
size_t i=0; i<coorddim; i++)
218 result[i] = 0.5 * (lower_[i] + upper_[i]);
233 if (dim == coorddim) {
234 for (
size_t i=0; i<coorddim; i++)
235 result[i] = (k & (1<<i)) ? upper_[i] : lower_[i];
239 unsigned int mask = 1;
241 for (
size_t i=0; i<coorddim; i++) {
243 result[i] = lower_[i];
245 result[i] = (k & mask) ? upper_[i] : lower_[i];
259 if (dim == coorddim) {
260 for (
size_t i=0; i<dim; i++)
261 vol *= upper_[i] - lower_[i];
263 }
else if (dim != 0) {
264 for (
size_t i=0; i<coorddim; i++)
266 vol *= upper_[i] - lower_[i];
281 for (
size_t i=0; i<dim; i++)
282 jacobianTransposed.diagonal()[i] = upper_[i] - lower_[i];
292 for (
size_t i=0; i<coorddim; i++)
294 jacobianTransposed[lc++][i] = upper_[i] - lower_[i];
300 for (
size_t i=0; i<dim; i++)
301 jacobianInverseTransposed.diagonal()[i] = 1.0 / (upper_[i] - lower_[i]);
311 for (
size_t i=0; i<coorddim; i++)
313 jacobianInverseTransposed[i][lc++] = 1.0 / (upper_[i] - lower_[i]);
316 Dune::FieldVector<ctype,coorddim> lower_;
318 Dune::FieldVector<ctype,coorddim> upper_;
320 std::bitset<coorddim> axes_;
const JacobianTransposed & jacobianTransposed(DUNE_UNUSED const LocalCoordinate &local) const
Jacobian transposed of the transformation from local to global coordinates.
Definition: axisalignedcubegeometry.hh:188
conditional< dim==coorddim, DiagonalMatrix< ctype, dim >, FieldMatrix< ctype, coorddim, dim > >::type JacobianInverseTransposed
Return type of jacobianInverseTransposed.
Definition: axisalignedcubegeometry.hh:88
A unique label for each type of element that can occur in a grid.
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
ctype volume() const
Return the element volume.
Definition: axisalignedcubegeometry.hh:256
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:129
const JacobianInverseTransposed & jacobianInverseTransposed(DUNE_UNUSED const LocalCoordinate &local) const
Jacobian transposed of the transformation from local to global coordinates.
Definition: axisalignedcubegeometry.hh:195
Definition: axisalignedcubegeometry.hh:59
GlobalCoordinate corner(int k) const
Return world coordinates of the k-th corner of the element.
Definition: axisalignedcubegeometry.hh:230
AxisAlignedCubeGeometry & operator=(const AxisAlignedCubeGeometry &other)
Assignment operator.
Definition: axisalignedcubegeometry.hh:136
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:210
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
std::size_t corners() const
Return the number of corners of the element.
Definition: axisalignedcubegeometry.hh:224
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:153
Cube element in any nonnegative dimension.
Definition: type.hh:31
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:24
GeometryType type() const
Type of the cube. Here: a hypercube of the correct dimension.
Definition: axisalignedcubegeometry.hh:147
LocalCoordinate local(const GlobalCoordinate &global) const
Map a point in global (world) coordinates to element coordinates.
Definition: axisalignedcubegeometry.hh:172
bool affine() const
Return if the element is affine. Here: yes.
Definition: axisalignedcubegeometry.hh:272
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:204
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