Dune Core Modules (2.7.1)

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 
13 #include <dune/common/fvector.hh>
14 #include <dune/common/fmatrix.hh>
16 #include <dune/common/unused.hh>
17 
18 #include <dune/geometry/referenceelements.hh>
19 #include <dune/geometry/type.hh>
20 
21 
22 namespace 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 
70  typedef ctype Volume;
71 
78  typedef typename std::conditional<dim==coorddim,
81 
88  typedef typename std::conditional<dim==coorddim,
91 
98  : lower_(lower),
99  upper_(upper),
100  axes_()
101  {
102  // all 'true', but is never actually used
103  axes_ = (1<<coorddim)-1;
104  }
105 
115  const std::bitset<coorddim>& axes)
116  : lower_(lower),
117  upper_(upper),
118  axes_(axes)
119  {
120  assert(axes.count()==dim);
121  for (size_t i=0; i<coorddim; i++)
122  if (not axes_[i])
123  upper_[i] = lower_[i];
124  }
125 
131  : lower_(lower)
132  {}
133 
136  {
137  lower_ = other.lower_;
138  upper_ = other.upper_;
139  axes_ = other.axes_;
140  return *this;
141  }
142 
145  {
146  return GeometryTypes::cube(dim);
147  }
148 
151  {
152  GlobalCoordinate result;
153  if (dim == coorddim) { // fast case
154  for (size_t i=0; i<coorddim; i++)
155  result[i] = lower_[i] + local[i]*(upper_[i] - lower_[i]);
156  } else if (dim == 0) { // a vertex -- the other fast case
157  result = lower_; // hope for named-return-type-optimization
158  } else { // slow case
159  size_t lc=0;
160  for (size_t i=0; i<coorddim; i++)
161  result[i] = (axes_[i])
162  ? lower_[i] + local[lc++]*(upper_[i] - lower_[i])
163  : lower_[i];
164  }
165  return result;
166  }
167 
170  {
171  LocalCoordinate result;
172  if (dim == coorddim) { // fast case
173  for (size_t i=0; i<dim; i++)
174  result[i] = (global[i] - lower_[i]) / (upper_[i] - lower_[i]);
175  } else if (dim != 0) { // slow case
176  size_t lc=0;
177  for (size_t i=0; i<coorddim; i++)
178  if (axes_[i])
179  result[lc++] = (global[i] - lower_[i]) / (upper_[i] - lower_[i]);
180  }
181  return result;
182  }
183 
186  {
187  JacobianTransposed result;
188 
189  // Actually compute the result. Uses different methods depending
190  // on what kind of matrix JacobianTransposed is.
191  jacobianTransposed(result);
192 
193  return result;
194  }
195 
198  {
200 
201  // Actually compute the result. Uses different methods depending
202  // on what kind of matrix JacobianTransposed is.
204 
205  return result;
206  }
207 
212  {
213  return volume();
214  }
215 
218  {
219  GlobalCoordinate result;
220  if (dim==0)
221  result = lower_;
222  else {
223  // Since lower_==upper_ for unused coordinates, this always does the right thing
224  for (size_t i=0; i<coorddim; i++)
225  result[i] = CoordType(0.5) * (lower_[i] + upper_[i]);
226  }
227  return result;
228  }
229 
231  int corners() const
232  {
233  return 1<<dim;
234  }
235 
238  {
239  GlobalCoordinate result;
240  if (dim == coorddim) { // fast case
241  for (size_t i=0; i<coorddim; i++)
242  result[i] = (k & (1<<i)) ? upper_[i] : lower_[i];
243  } else if (dim == 0) { // vertex
244  result = lower_; // rely on named return-type optimization
245  } else { // slow case
246  unsigned int mask = 1;
247 
248  for (size_t i=0; i<coorddim; i++) {
249  if (not axes_[i])
250  result[i] = lower_[i];
251  else {
252  result[i] = (k & mask) ? upper_[i] : lower_[i];
253  mask = (mask<<1);
254  }
255  }
256  }
257 
258 
259  return result;
260  }
261 
263  Volume volume() const
264  {
265  ctype vol = 1;
266  if (dim == coorddim) { // fast case
267  for (size_t i=0; i<dim; i++)
268  vol *= upper_[i] - lower_[i];
269  // do nothing if dim == 0
270  } else if (dim != 0) { // slow case
271  for (size_t i=0; i<coorddim; i++)
272  if (axes_[i])
273  vol *= upper_[i] - lower_[i];
274  }
275  return vol;
276  }
277 
279  bool affine() const
280  {
281  return true;
282  }
283 
284  friend Dune::Transitional::ReferenceElement< ctype, Dim<dim> > referenceElement ( const AxisAlignedCubeGeometry &geometry )
285  {
287  }
288 
289  private:
290  // jacobianTransposed: fast case --> diagonal matrix
291  void jacobianTransposed ( DiagonalMatrix<ctype,dim> &jacobianTransposed ) const
292  {
293  for (size_t i=0; i<dim; i++)
294  jacobianTransposed.diagonal()[i] = upper_[i] - lower_[i];
295  }
296 
297  // jacobianTransposed: slow case --> dense matrix
298  void jacobianTransposed ( FieldMatrix<ctype,dim,coorddim> &jacobianTransposed ) const
299  {
300  if (dim==0)
301  return;
302 
303  size_t lc = 0;
304  for (size_t i=0; i<coorddim; i++)
305  if (axes_[i])
306  jacobianTransposed[lc++][i] = upper_[i] - lower_[i];
307  }
308 
309  // jacobianInverseTransposed: fast case --> diagonal matrix
310  void jacobianInverseTransposed ( DiagonalMatrix<ctype,dim> &jacobianInverseTransposed ) const
311  {
312  for (size_t i=0; i<dim; i++)
313  jacobianInverseTransposed.diagonal()[i] = CoordType(1.0) / (upper_[i] - lower_[i]);
314  }
315 
316  // jacobianInverseTransposed: slow case --> dense matrix
317  void jacobianInverseTransposed ( FieldMatrix<ctype,coorddim,dim> &jacobianInverseTransposed ) const
318  {
319  if (dim==0)
320  return;
321 
322  size_t lc = 0;
323  for (size_t i=0; i<coorddim; i++)
324  if (axes_[i])
325  jacobianInverseTransposed[i][lc++] = CoordType(1.0) / (upper_[i] - lower_[i]);
326  }
327 
329 
331 
332  std::bitset<coorddim> axes_;
333  };
334 
335 } // namespace Dune
336 #endif
A geometry implementation for axis-aligned hypercubes.
Definition: axisalignedcubegeometry.hh:49
Volume volume() const
Return the element volume.
Definition: axisalignedcubegeometry.hh:263
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:113
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:96
AxisAlignedCubeGeometry(const Dune::FieldVector< ctype, coorddim > lower)
Constructor from a single point only.
Definition: axisalignedcubegeometry.hh:130
GlobalCoordinate corner(int k) const
Return world coordinates of the k-th corner of the element.
Definition: axisalignedcubegeometry.hh:237
JacobianTransposed jacobianTransposed(DUNE_UNUSED const LocalCoordinate &local) const
Jacobian transposed of the transformation from local to global coordinates.
Definition: axisalignedcubegeometry.hh:185
ctype Volume
Type used for volume.
Definition: axisalignedcubegeometry.hh:70
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:135
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:197
LocalCoordinate local(const GlobalCoordinate &global) const
Map a point in global (world) coordinates to element coordinates.
Definition: axisalignedcubegeometry.hh:169
CoordType ctype
Type used for single coordinate coefficients.
Definition: axisalignedcubegeometry.hh:61
std::conditional< dim==coorddim, DiagonalMatrix< ctype, dim >, FieldMatrix< ctype, dim, coorddim > >::type JacobianTransposed
Return type of jacobianTransposed.
Definition: axisalignedcubegeometry.hh:80
GeometryType type() const
Type of the cube. Here: a hypercube of the correct dimension.
Definition: axisalignedcubegeometry.hh:144
int corners() const
Return the number of corners of the element.
Definition: axisalignedcubegeometry.hh:231
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:211
std::conditional< dim==coorddim, DiagonalMatrix< ctype, dim >, FieldMatrix< ctype, coorddim, dim > >::type JacobianInverseTransposed
Return type of jacobianInverseTransposed.
Definition: axisalignedcubegeometry.hh:90
GlobalCoordinate center() const
Return center of mass of the element.
Definition: axisalignedcubegeometry.hh:217
bool affine() const
Return if the element is affine. Here: yes.
Definition: axisalignedcubegeometry.hh:279
GlobalCoordinate global(const LocalCoordinate &local) const
Map a point in local (element) coordinates to world coordinates.
Definition: axisalignedcubegeometry.hh:150
A diagonal matrix of static size.
Definition: diagonalmatrix.hh:52
A dense n x m matrix.
Definition: fmatrix.hh:69
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:279
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:775
Mask< V > mask(ADLTag< 0, std::is_same< V, Mask< V > >::value >, const V &v)
implements Simd::mask()
Definition: defaults.hh:151
Dune namespace.
Definition: alignedallocator.hh:14
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.80.0 (May 16, 22:29, 2024)