Dune Core Modules (2.4.2)

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/type.hh>
19 
20 
21 
22 
23 namespace Dune {
24 
48  template <class CoordType, unsigned int dim, unsigned int coorddim>
50  {
51 
52 
53  public:
54 
56  enum {mydimension = dim};
57 
59  enum {coorddimension = coorddim};
60 
62  typedef CoordType ctype;
63 
66 
69 
76  typedef typename conditional<dim==coorddim,
79 
86  typedef typename conditional<dim==coorddim,
89 
96  : lower_(lower),
97  upper_(upper),
98  axes_()
99  {
100  // all 'true', but is never actually used
101  axes_ = (1<<coorddim)-1;
102  }
103 
113  const std::bitset<coorddim>& axes)
114  : lower_(lower),
115  upper_(upper),
116  axes_(axes)
117  {
118  assert(axes.count()==dim);
119  for (size_t i=0; i<coorddim; i++)
120  if (not axes_[i])
121  upper_[i] = lower_[i];
122  }
123 
129  : lower_(lower)
130  {}
131 
134  {
135  lower_ = other.lower_;
136  upper_ = other.upper_;
137  axes_ = other.axes_;
138  return *this;
139  }
140 
143  {
144  return GeometryType(GeometryType::cube,dim);
145  }
146 
149  {
150  GlobalCoordinate result;
151  if (dim == coorddim) { // fast case
152  for (size_t i=0; i<coorddim; i++)
153  result[i] = lower_[i] + local[i]*(upper_[i] - lower_[i]);
154  } if (dim == 0) { // a vertex -- the other fast case
155  result = lower_; // hope for named-return-type-optimization
156  } else { // slow case
157  size_t lc=0;
158  for (size_t i=0; i<coorddim; i++)
159  result[i] = (axes_[i])
160  ? lower_[i] + local[lc++]*(upper_[i] - lower_[i])
161  : lower_[i];
162  }
163  return result;
164  }
165 
168  {
169  LocalCoordinate result;
170  if (dim == coorddim) { // fast case
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) { // slow case
174  size_t lc=0;
175  for (size_t i=0; i<coorddim; i++)
176  if (axes_[i])
177  result[lc++] = (global[i] - lower_[i]) / (upper_[i] - lower_[i]);
178  }
179  return result;
180  }
181 
184  {
185  JacobianTransposed result;
186 
187  // Actually compute the result. Computes different methods depending
188  // on what kind of matrix JacobianTransposed is.
189  jacobianTransposed(result);
190 
191  return result;
192  }
193 
196  {
198 
199  // Actually compute the result. Computes different methods depending
200  // on what kind of matrix JacobianTransposed is.
202 
203  return result;
204  }
205 
209  ctype integrationElement(DUNE_UNUSED const LocalCoordinate& local) const
210  {
211  return volume();
212  }
213 
216  {
217  GlobalCoordinate result;
218  if (dim==0)
219  result = lower_;
220  else {
221  // Since lower_==upper_ for unused coordinates, this always does the right thing
222  for (size_t i=0; i<coorddim; i++)
223  result[i] = 0.5 * (lower_[i] + upper_[i]);
224  }
225  return result;
226  }
227 
229  int corners() const
230  {
231  return 1<<dim;
232  }
233 
236  {
237  GlobalCoordinate result;
238  if (dim == coorddim) { // fast case
239  for (size_t i=0; i<coorddim; i++)
240  result[i] = (k & (1<<i)) ? upper_[i] : lower_[i];
241  } if (dim == 0) { // vertex
242  result = lower_; // rely on named return-type optimization
243  } else { // slow case
244  unsigned int mask = 1;
245 
246  for (size_t i=0; i<coorddim; i++) {
247  if (not axes_[i])
248  result[i] = lower_[i];
249  else {
250  result[i] = (k & mask) ? upper_[i] : lower_[i];
251  mask = (mask<<1);
252  }
253  }
254  }
255 
256 
257  return result;
258  }
259 
261  ctype volume() const
262  {
263  ctype vol = 1;
264  if (dim == coorddim) { // fast case
265  for (size_t i=0; i<dim; i++)
266  vol *= upper_[i] - lower_[i];
267  // do nothing if dim == 0
268  } else if (dim != 0) { // slow case
269  for (size_t i=0; i<coorddim; i++)
270  if (axes_[i])
271  vol *= upper_[i] - lower_[i];
272  }
273  return vol;
274  }
275 
277  bool affine() const
278  {
279  return true;
280  }
281 
282  private:
283  // jacobianTransposed: fast case --> diagonal matrix
285  {
286  for (size_t i=0; i<dim; i++)
287  jacobianTransposed.diagonal()[i] = upper_[i] - lower_[i];
288  }
289 
290  // jacobianTransposed: slow case --> dense matrix
291  void jacobianTransposed ( FieldMatrix<ctype,dim,coorddim> &jacobianTransposed ) const
292  {
293  if (dim==0)
294  return;
295 
296  size_t lc = 0;
297  for (size_t i=0; i<coorddim; i++)
298  if (axes_[i])
299  jacobianTransposed[lc++][i] = upper_[i] - lower_[i];
300  }
301 
302  // jacobianInverseTransposed: fast case --> diagonal matrix
303  void jacobianInverseTransposed ( DiagonalMatrix<ctype,dim> &jacobianInverseTransposed ) const
304  {
305  for (size_t i=0; i<dim; i++)
306  jacobianInverseTransposed.diagonal()[i] = 1.0 / (upper_[i] - lower_[i]);
307  }
308 
309  // jacobianInverseTransposed: slow case --> dense matrix
310  void jacobianInverseTransposed ( FieldMatrix<ctype,coorddim,dim> &jacobianInverseTransposed ) const
311  {
312  if (dim==0)
313  return;
314 
315  size_t lc = 0;
316  for (size_t i=0; i<coorddim; i++)
317  if (axes_[i])
318  jacobianInverseTransposed[i][lc++] = 1.0 / (upper_[i] - lower_[i]);
319  }
320 
322 
324 
325  std::bitset<coorddim> axes_;
326  };
327 
328 } // namespace Dune
329 #endif
A geometry implementation for axis-aligned hypercubes.
Definition: axisalignedcubegeometry.hh:50
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
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
AxisAlignedCubeGeometry(const Dune::FieldVector< ctype, coorddim > lower)
Constructor from a single point only.
Definition: axisalignedcubegeometry.hh:128
conditional< dim==coorddim, DiagonalMatrix< ctype, dim >, FieldMatrix< ctype, coorddim, dim > >::type JacobianInverseTransposed
Return type of jacobianInverseTransposed.
Definition: axisalignedcubegeometry.hh:88
GlobalCoordinate corner(int k) const
Return world coordinates of the k-th corner of the element.
Definition: axisalignedcubegeometry.hh:235
JacobianTransposed jacobianTransposed(DUNE_UNUSED const LocalCoordinate &local) const
Jacobian transposed of the transformation from local to global coordinates.
Definition: axisalignedcubegeometry.hh:183
FieldVector< ctype, dim > LocalCoordinate
Type used for a vector of element coordinates.
Definition: axisalignedcubegeometry.hh:65
AxisAlignedCubeGeometry & operator=(const AxisAlignedCubeGeometry &other)
Assignment operator.
Definition: axisalignedcubegeometry.hh:133
FieldVector< ctype, coorddim > GlobalCoordinate
Type used for a vector of world coordinates.
Definition: axisalignedcubegeometry.hh:68
JacobianInverseTransposed jacobianInverseTransposed(DUNE_UNUSED const LocalCoordinate &local) const
Jacobian transposed of the transformation from local to global coordinates.
Definition: axisalignedcubegeometry.hh:195
LocalCoordinate local(const GlobalCoordinate &global) const
Map a point in global (world) coordinates to element coordinates.
Definition: axisalignedcubegeometry.hh:167
CoordType ctype
Type used for single coordinate coefficients.
Definition: axisalignedcubegeometry.hh:62
GeometryType type() const
Type of the cube. Here: a hypercube of the correct dimension.
Definition: axisalignedcubegeometry.hh:142
int corners() const
Return the number of corners of the element.
Definition: axisalignedcubegeometry.hh:229
conditional< dim==coorddim, DiagonalMatrix< ctype, dim >, FieldMatrix< ctype, dim, coorddim > >::type JacobianTransposed
Return type of jacobianTransposed.
Definition: axisalignedcubegeometry.hh:78
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
GlobalCoordinate center() const
Return center of mass of the element.
Definition: axisalignedcubegeometry.hh:215
ctype volume() const
Return the element volume.
Definition: axisalignedcubegeometry.hh:261
bool affine() const
Return if the element is affine. Here: yes.
Definition: axisalignedcubegeometry.hh:277
GlobalCoordinate global(const LocalCoordinate &local) const
Map a point in local (element) coordinates to world coordinates.
Definition: axisalignedcubegeometry.hh:148
A diagonal matrix of static size.
Definition: diagonalmatrix.hh:51
A dense n x m matrix.
Definition: fmatrix.hh:67
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
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.
Dune namespace.
Definition: alignment.hh:10
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)