Dune Core Modules (2.6.0)

axisalignedcubegeometry.hh
Go to the documentation of this file.
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3
4#ifndef DUNE_GEOMETRY_AXISALIGNED_CUBE_GEOMETRY_HH
5#define DUNE_GEOMETRY_AXISALIGNED_CUBE_GEOMETRY_HH
6
11#include <bitset>
12
16#include <dune/common/unused.hh>
17
18#include <dune/geometry/referenceelements.hh>
19#include <dune/geometry/type.hh>
20
21
22namespace Dune {
23
47 template <class CoordType, unsigned int dim, unsigned int coorddim>
49 {
50
51
52 public:
53
55 enum {mydimension = dim};
56
58 enum {coorddimension = coorddim};
59
61 typedef CoordType ctype;
62
65
68
75 typedef typename std::conditional<dim==coorddim,
78
85 typedef typename std::conditional<dim==coorddim,
88
95 : lower_(lower),
96 upper_(upper),
97 axes_()
98 {
99 // all 'true', but is never actually used
100 axes_ = (1<<coorddim)-1;
101 }
102
112 const std::bitset<coorddim>& axes)
113 : lower_(lower),
114 upper_(upper),
115 axes_(axes)
116 {
117 assert(axes.count()==dim);
118 for (size_t i=0; i<coorddim; i++)
119 if (not axes_[i])
120 upper_[i] = lower_[i];
121 }
122
128 : lower_(lower)
129 {}
130
133 {
134 lower_ = other.lower_;
135 upper_ = other.upper_;
136 axes_ = other.axes_;
137 return *this;
138 }
139
142 {
143 return GeometryTypes::cube(dim);
144 }
145
148 {
149 GlobalCoordinate result;
150 if (dim == coorddim) { // fast case
151 for (size_t i=0; i<coorddim; i++)
152 result[i] = lower_[i] + local[i]*(upper_[i] - lower_[i]);
153 } else if (dim == 0) { // a vertex -- the other fast case
154 result = lower_; // hope for named-return-type-optimization
155 } else { // slow case
156 size_t lc=0;
157 for (size_t i=0; i<coorddim; i++)
158 result[i] = (axes_[i])
159 ? lower_[i] + local[lc++]*(upper_[i] - lower_[i])
160 : lower_[i];
161 }
162 return result;
163 }
164
167 {
168 LocalCoordinate result;
169 if (dim == coorddim) { // fast case
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) { // slow case
173 size_t lc=0;
174 for (size_t i=0; i<coorddim; i++)
175 if (axes_[i])
176 result[lc++] = (global[i] - lower_[i]) / (upper_[i] - lower_[i]);
177 }
178 return result;
179 }
180
183 {
184 JacobianTransposed result;
185
186 // Actually compute the result. Uses different methods depending
187 // on what kind of matrix JacobianTransposed is.
188 jacobianTransposed(result);
189
190 return result;
191 }
192
195 {
197
198 // Actually compute the result. Uses different methods depending
199 // on what kind of matrix JacobianTransposed is.
201
202 return result;
203 }
204
209 {
210 return volume();
211 }
212
215 {
216 GlobalCoordinate result;
217 if (dim==0)
218 result = lower_;
219 else {
220 // Since lower_==upper_ for unused coordinates, this always does the right thing
221 for (size_t i=0; i<coorddim; i++)
222 result[i] = 0.5 * (lower_[i] + upper_[i]);
223 }
224 return result;
225 }
226
228 int corners() const
229 {
230 return 1<<dim;
231 }
232
235 {
236 GlobalCoordinate result;
237 if (dim == coorddim) { // fast case
238 for (size_t i=0; i<coorddim; i++)
239 result[i] = (k & (1<<i)) ? upper_[i] : lower_[i];
240 } else if (dim == 0) { // vertex
241 result = lower_; // rely on named return-type optimization
242 } else { // slow case
243 unsigned int mask = 1;
244
245 for (size_t i=0; i<coorddim; i++) {
246 if (not axes_[i])
247 result[i] = lower_[i];
248 else {
249 result[i] = (k & mask) ? upper_[i] : lower_[i];
250 mask = (mask<<1);
251 }
252 }
253 }
254
255
256 return result;
257 }
258
260 ctype volume() const
261 {
262 ctype vol = 1;
263 if (dim == coorddim) { // fast case
264 for (size_t i=0; i<dim; i++)
265 vol *= upper_[i] - lower_[i];
266 // do nothing if dim == 0
267 } else if (dim != 0) { // slow case
268 for (size_t i=0; i<coorddim; i++)
269 if (axes_[i])
270 vol *= upper_[i] - lower_[i];
271 }
272 return vol;
273 }
274
276 bool affine() const
277 {
278 return true;
279 }
280
281 friend Dune::Transitional::ReferenceElement< ctype, Dim<dim> > referenceElement ( const AxisAlignedCubeGeometry &geometry )
282 {
284 }
285
286 private:
287 // jacobianTransposed: fast case --> diagonal matrix
288 void jacobianTransposed ( DiagonalMatrix<ctype,dim> &jacobianTransposed ) const
289 {
290 for (size_t i=0; i<dim; i++)
291 jacobianTransposed.diagonal()[i] = upper_[i] - lower_[i];
292 }
293
294 // jacobianTransposed: slow case --> dense matrix
295 void jacobianTransposed ( FieldMatrix<ctype,dim,coorddim> &jacobianTransposed ) const
296 {
297 if (dim==0)
298 return;
299
300 size_t lc = 0;
301 for (size_t i=0; i<coorddim; i++)
302 if (axes_[i])
303 jacobianTransposed[lc++][i] = upper_[i] - lower_[i];
304 }
305
306 // jacobianInverseTransposed: fast case --> diagonal matrix
307 void jacobianInverseTransposed ( DiagonalMatrix<ctype,dim> &jacobianInverseTransposed ) const
308 {
309 for (size_t i=0; i<dim; i++)
310 jacobianInverseTransposed.diagonal()[i] = 1.0 / (upper_[i] - lower_[i]);
311 }
312
313 // jacobianInverseTransposed: slow case --> dense matrix
314 void jacobianInverseTransposed ( FieldMatrix<ctype,coorddim,dim> &jacobianInverseTransposed ) const
315 {
316 if (dim==0)
317 return;
318
319 size_t lc = 0;
320 for (size_t i=0; i<coorddim; i++)
321 if (axes_[i])
322 jacobianInverseTransposed[i][lc++] = 1.0 / (upper_[i] - lower_[i]);
323 }
324
326
328
329 std::bitset<coorddim> axes_;
330 };
331
332} // namespace Dune
333#endif
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
std::conditional< dim==coorddim, DiagonalMatrix< ctype, dim >, FieldMatrix< ctype, dim, coorddim > >::type JacobianTransposed
Return type of jacobianTransposed.
Definition: axisalignedcubegeometry.hh:77
AxisAlignedCubeGeometry(const Dune::FieldVector< ctype, coorddim > lower)
Constructor from a single point only.
Definition: axisalignedcubegeometry.hh:127
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
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
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
std::conditional< dim==coorddim, DiagonalMatrix< ctype, dim >, FieldMatrix< ctype, coorddim, dim > >::type JacobianInverseTransposed
Return type of jacobianInverseTransposed.
Definition: axisalignedcubegeometry.hh:87
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
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
A diagonal matrix of static size.
Definition: diagonalmatrix.hh:52
A dense n x m matrix.
Definition: fmatrix.hh:68
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:277
This file implements a quadratic diagonal matrix of fixed size.
Implements a matrix constructed from a given type representing a field and compile-time given number ...
Implements a vector constructed from a given type representing a field and a compile-time given size.
#define DUNE_UNUSED
A macro for marking variables that the compiler mistakenly flags as unused, which sometimes happens d...
Definition: unused.hh:16
unspecified-type ReferenceElement
Returns the type of reference element for the argument types T...
Definition: referenceelements.hh:415
constexpr GeometryType cube(unsigned int dim)
Returns a GeometryType representing a hypercube of dimension dim.
Definition: type.hh:705
Dune namespace.
Definition: alignedallocator.hh:10
static const ReferenceElement & cube()
get hypercube reference elements
Definition: referenceelements.hh:208
A unique label for each type of element that can occur in a grid.
Definition of the DUNE_UNUSED macro for the case that config.h is not available.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 13, 23:29, 2024)