dune-geometry  2.3beta2
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>
15 #include <dune/common/diagonalmatrix.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 
65  typedef FieldVector<ctype,dim> LocalCoordinate;
66 
68  typedef FieldVector<ctype,coorddim> GlobalCoordinate;
69 
76  typedef typename conditional<dim==coorddim,
77  DiagonalMatrix<ctype,dim>,
78  FieldMatrix<ctype,dim,coorddim> >::type JacobianTransposed;
79 
86  typedef typename conditional<dim==coorddim,
87  DiagonalMatrix<ctype,dim>,
88  FieldMatrix<ctype,coorddim,dim> >::type JacobianInverseTransposed;
89 
94  AxisAlignedCubeGeometry(const Dune::FieldVector<ctype,coorddim> lower,
95  const Dune::FieldVector<ctype,coorddim> upper)
96  : lower_(lower),
97  upper_(upper),
98  axes_((1<<coorddim)-1), // all 'true', but is never actually used
99  jacobianTransposed_(0),
100  jacobianInverseTransposed_(0)
101  {}
102 
110  AxisAlignedCubeGeometry(const Dune::FieldVector<ctype,coorddim> lower,
111  const Dune::FieldVector<ctype,coorddim> upper,
112  const std::bitset<coorddim>& axes)
113  : lower_(lower),
114  upper_(upper),
115  axes_(axes),
116  jacobianTransposed_(0),
117  jacobianInverseTransposed_(0)
118  {
119  assert(axes.count()==dim);
120  for (size_t i=0; i<coorddim; i++)
121  if (not axes_[i])
122  upper_[i] = lower_[i];
123  }
124 
129  AxisAlignedCubeGeometry(const Dune::FieldVector<ctype,coorddim> lower)
130  : lower_(lower),
131  jacobianTransposed_(0),
132  jacobianInverseTransposed_(0)
133  {}
134 
137  {
138  lower_ = other.lower_;
139  upper_ = other.upper_;
140  axes_ = other.axes_;
141  jacobianTransposed_ = other.jacobianTransposed_;
142  jacobianInverseTransposed_ = other.jacobianInverseTransposed_;
143  return *this;
144  }
145 
148  {
149  return GeometryType(GeometryType::cube,dim);
150  }
151 
154  {
155  GlobalCoordinate result;
156  if (dim == coorddim) { // fast case
157  for (size_t i=0; i<coorddim; i++)
158  result[i] = lower_[i] + local[i]*(upper_[i] - lower_[i]);
159  } if (dim == 0) { // a vertex -- the other fast case
160  result = lower_; // hope for named-return-type-optimization
161  } else { // slow case
162  size_t lc=0;
163  for (size_t i=0; i<coorddim; i++)
164  result[i] = (axes_[i])
165  ? lower_[i] + local[lc++]*(upper_[i] - lower_[i])
166  : lower_[i];
167  }
168  return result;
169  }
170 
173  {
174  LocalCoordinate result;
175  if (dim == coorddim) { // fast case
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) { // slow case
179  size_t lc=0;
180  for (size_t i=0; i<coorddim; i++)
181  if (axes_[i])
182  result[lc++] = (global[i] - lower_[i]) / (upper_[i] - lower_[i]);
183  }
184  return result;
185  }
186 
189  {
190  jacobianTransposed( jacobianTransposed_ );
191  return jacobianTransposed_;
192  }
193 
196  {
197  jacobianInverseTransposed( jacobianInverseTransposed_ );
198  return jacobianInverseTransposed_;
199  }
200 
204  ctype integrationElement(DUNE_UNUSED const LocalCoordinate& local) const
205  {
206  return volume();
207  }
208 
211  {
212  GlobalCoordinate result;
213  if (dim==0)
214  result = lower_;
215  else {
216  // Since lower_==upper_ for unused coordinates, this always does the right thing
217  for (size_t i=0; i<coorddim; i++)
218  result[i] = 0.5 * (lower_[i] + upper_[i]);
219  }
220  return result;
221  }
222 
224  std::size_t corners() const
225  {
226  return 1<<dim;
227  }
228 
231  {
232  GlobalCoordinate result;
233  if (dim == coorddim) { // fast case
234  for (size_t i=0; i<coorddim; i++)
235  result[i] = (k & (1<<i)) ? upper_[i] : lower_[i];
236  } if (dim == 0) { // vertex
237  result = lower_; // rely on named return-type optimization
238  } else { // slow case
239  unsigned int mask = 1;
240 
241  for (size_t i=0; i<coorddim; i++) {
242  if (not axes_[i])
243  result[i] = lower_[i];
244  else {
245  result[i] = (k & mask) ? upper_[i] : lower_[i];
246  mask = (mask<<1);
247  }
248  }
249  }
250 
251 
252  return result;
253  }
254 
256  ctype volume() const
257  {
258  ctype vol = 1;
259  if (dim == coorddim) { // fast case
260  for (size_t i=0; i<dim; i++)
261  vol *= upper_[i] - lower_[i];
262  // do nothing if dim == 0
263  } else if (dim != 0) { // slow case
264  for (size_t i=0; i<coorddim; i++)
265  if (axes_[i])
266  vol *= upper_[i] - lower_[i];
267  }
268  return vol;
269  }
270 
272  bool affine() const
273  {
274  return true;
275  }
276 
277  private:
278  // jacobianTransposed: fast case --> diagonal matrix
279  void jacobianTransposed ( DiagonalMatrix<ctype,dim> &jacobianTransposed ) const
280  {
281  for (size_t i=0; i<dim; i++)
282  jacobianTransposed.diagonal()[i] = upper_[i] - lower_[i];
283  }
284 
285  // jacobianTransposed: slow case --> dense matrix
286  void jacobianTransposed ( FieldMatrix<ctype,dim,coorddim> &jacobianTransposed ) const
287  {
288  if (dim==0)
289  return;
290 
291  size_t lc = 0;
292  for (size_t i=0; i<coorddim; i++)
293  if (axes_[i])
294  jacobianTransposed[lc++][i] = upper_[i] - lower_[i];
295  }
296 
297  // jacobianInverseTransposed: fast case --> diagonal matrix
298  void jacobianInverseTransposed ( DiagonalMatrix<ctype,dim> &jacobianInverseTransposed ) const
299  {
300  for (size_t i=0; i<dim; i++)
301  jacobianInverseTransposed.diagonal()[i] = 1.0 / (upper_[i] - lower_[i]);
302  }
303 
304  // jacobianInverseTransposed: slow case --> dense matrix
305  void jacobianInverseTransposed ( FieldMatrix<ctype,coorddim,dim> &jacobianInverseTransposed ) const
306  {
307  if (dim==0)
308  return;
309 
310  size_t lc = 0;
311  for (size_t i=0; i<coorddim; i++)
312  if (axes_[i])
313  jacobianInverseTransposed[i][lc++] = 1.0 / (upper_[i] - lower_[i]);
314  }
315 
316  Dune::FieldVector<ctype,coorddim> lower_;
317 
318  Dune::FieldVector<ctype,coorddim> upper_;
319 
320  std::bitset<coorddim> axes_;
321 
322  // Storage so method jacobianTransposed can return a const reference
323  mutable JacobianTransposed jacobianTransposed_;
324 
325  // Storage so method jacobianInverseTransposed can return a const reference
326  mutable JacobianInverseTransposed jacobianInverseTransposed_;
327 
328  };
329 
330 } // namespace Dune
331 #endif