Dune Core Modules (2.9.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// SPDX-FileCopyrightInfo: Copyright (C) DUNE Project contributors, see file LICENSE.md in module root
4// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
5
6#ifndef DUNE_GEOMETRY_AXISALIGNED_CUBE_GEOMETRY_HH
7#define DUNE_GEOMETRY_AXISALIGNED_CUBE_GEOMETRY_HH
8
13#include <bitset>
14
18
19#include <dune/geometry/referenceelements.hh>
20#include <dune/geometry/type.hh>
21
22
23namespace Dune {
24
48 template <class CoordType, unsigned int dim, unsigned int coorddim>
50 {
51
52
53 public:
54
56 constexpr static int mydimension = dim;
57
59 constexpr static int coorddimension = coorddim;
60
62 typedef CoordType ctype;
63
66
69
71 typedef ctype Volume;
72
79 typedef typename std::conditional<dim==coorddim,
82
89 typedef typename std::conditional<dim==coorddim,
92
100 using Jacobian = std::conditional_t<dim==coorddim, DiagonalMatrix<ctype,dim>, FieldMatrix<ctype,coorddim,dim> >;
101
109 using JacobianInverse = std::conditional_t<dim==coorddim, DiagonalMatrix<ctype,dim>, FieldMatrix<ctype,dim,coorddim> >;
110
117 : lower_(lower),
118 upper_(upper),
119 axes_()
120 {
121 static_assert(dim==coorddim, "Use this constructor only if dim==coorddim!");
122 // all 'true', but is never actually used
123 axes_ = (1<<coorddim)-1;
124 }
125
135 const std::bitset<coorddim>& axes)
136 : lower_(lower),
137 upper_(upper),
138 axes_(axes)
139 {
140 assert(axes.count()==dim);
141 for (size_t i=0; i<coorddim; i++)
142 if (not axes_[i])
143 upper_[i] = lower_[i];
144 }
145
151 : lower_(lower)
152 {}
153
156 {
157 return GeometryTypes::cube(dim);
158 }
159
162 {
163 GlobalCoordinate result;
164 if (dim == coorddim) { // fast case
165 for (size_t i=0; i<coorddim; i++)
166 result[i] = lower_[i] + local[i]*(upper_[i] - lower_[i]);
167 } else if (dim == 0) { // a vertex -- the other fast case
168 result = lower_; // hope for named-return-type-optimization
169 } else { // slow case
170 size_t lc=0;
171 for (size_t i=0; i<coorddim; i++)
172 result[i] = (axes_[i])
173 ? lower_[i] + local[lc++]*(upper_[i] - lower_[i])
174 : lower_[i];
175 }
176 return result;
177 }
178
181 {
182 LocalCoordinate result;
183 if (dim == coorddim) { // fast case
184 for (size_t i=0; i<dim; i++)
185 result[i] = (global[i] - lower_[i]) / (upper_[i] - lower_[i]);
186 } else if (dim != 0) { // slow case
187 size_t lc=0;
188 for (size_t i=0; i<coorddim; i++)
189 if (axes_[i])
190 result[lc++] = (global[i] - lower_[i]) / (upper_[i] - lower_[i]);
191 }
192 return result;
193 }
194
197 {
198 JacobianTransposed result;
199
200 // Actually compute the result. Uses different methods depending
201 // on what kind of matrix JacobianTransposed is.
202 jacobianTransposed(result);
203
204 return result;
205 }
206
209 {
211
212 // Actually compute the result. Uses different methods depending
213 // on what kind of matrix JacobianTransposed is.
215
216 return result;
217 }
218
220 Jacobian jacobian([[maybe_unused]] const LocalCoordinate& local) const
221 {
222 return jacobianTransposed(local).transposed();
223 }
224
227 {
228 return jacobianInverseTransposed(local).transposed();
229 }
230
234 Volume integrationElement([[maybe_unused]] const LocalCoordinate& local) const
235 {
236 return volume();
237 }
238
241 {
242 GlobalCoordinate result;
243 if (dim==0)
244 result = lower_;
245 else {
246 // Since lower_==upper_ for unused coordinates, this always does the right thing
247 for (size_t i=0; i<coorddim; i++)
248 result[i] = CoordType(0.5) * (lower_[i] + upper_[i]);
249 }
250 return result;
251 }
252
254 int corners() const
255 {
256 return 1<<dim;
257 }
258
261 {
262 GlobalCoordinate result;
263 if (dim == coorddim) { // fast case
264 for (size_t i=0; i<coorddim; i++)
265 result[i] = (k & (1<<i)) ? upper_[i] : lower_[i];
266 } else if (dim == 0) { // vertex
267 result = lower_; // rely on named return-type optimization
268 } else { // slow case
269 unsigned int mask = 1;
270
271 for (size_t i=0; i<coorddim; i++) {
272 if (not axes_[i])
273 result[i] = lower_[i];
274 else {
275 result[i] = (k & mask) ? upper_[i] : lower_[i];
276 mask = (mask<<1);
277 }
278 }
279 }
280
281
282 return result;
283 }
284
287 {
288 ctype vol = 1;
289 if (dim == coorddim) { // fast case
290 for (size_t i=0; i<dim; i++)
291 vol *= upper_[i] - lower_[i];
292 // do nothing if dim == 0
293 } else if (dim != 0) { // slow case
294 for (size_t i=0; i<coorddim; i++)
295 if (axes_[i])
296 vol *= upper_[i] - lower_[i];
297 }
298 return vol;
299 }
300
302 bool affine() const
303 {
304 return true;
305 }
306
307 friend Dune::Transitional::ReferenceElement< ctype, Dim<dim> > referenceElement ( const AxisAlignedCubeGeometry & /* geometry */ )
308 {
310 }
311
312 private:
313 // jacobianTransposed: fast case --> diagonal matrix
314 void jacobianTransposed ( DiagonalMatrix<ctype,dim> &jacobianTransposed ) const
315 {
316 for (size_t i=0; i<dim; i++)
317 jacobianTransposed.diagonal()[i] = upper_[i] - lower_[i];
318 }
319
320 // jacobianTransposed: slow case --> dense matrix
321 void jacobianTransposed ( FieldMatrix<ctype,dim,coorddim> &jacobianTransposed ) const
322 {
323 if (dim==0)
324 return;
325
326 size_t lc = 0;
327 for (size_t i=0; i<coorddim; i++)
328 if (axes_[i])
329 jacobianTransposed[lc++][i] = upper_[i] - lower_[i];
330 }
331
332 // jacobianInverseTransposed: fast case --> diagonal matrix
333 void jacobianInverseTransposed ( DiagonalMatrix<ctype,dim> &jacobianInverseTransposed ) const
334 {
335 for (size_t i=0; i<dim; i++)
336 jacobianInverseTransposed.diagonal()[i] = CoordType(1.0) / (upper_[i] - lower_[i]);
337 }
338
339 // jacobianInverseTransposed: slow case --> dense matrix
340 void jacobianInverseTransposed ( FieldMatrix<ctype,coorddim,dim> &jacobianInverseTransposed ) const
341 {
342 if (dim==0)
343 return;
344
345 size_t lc = 0;
346 for (size_t i=0; i<coorddim; i++)
347 if (axes_[i])
348 jacobianInverseTransposed[i][lc++] = CoordType(1.0) / (upper_[i] - lower_[i]);
349 }
350
352
354
355 std::bitset<coorddim> axes_;
356 };
357
358} // namespace Dune
359#endif
A geometry implementation for axis-aligned hypercubes.
Definition: axisalignedcubegeometry.hh:50
Volume volume() const
Return the element volume.
Definition: axisalignedcubegeometry.hh:286
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:133
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:115
JacobianInverse jacobianInverse(const LocalCoordinate &local) const
Inverse Jacobian of the transformation from local to global coordinates.
Definition: axisalignedcubegeometry.hh:226
JacobianInverseTransposed jacobianInverseTransposed(const LocalCoordinate &local) const
Inverse Jacobian transposed of the transformation from local to global coordinates.
Definition: axisalignedcubegeometry.hh:208
std::conditional< dim==coorddim, DiagonalMatrix< ctype, dim >, FieldMatrix< ctype, dim, coorddim > >::type JacobianTransposed
Return type of jacobianTransposed.
Definition: axisalignedcubegeometry.hh:81
AxisAlignedCubeGeometry(const Dune::FieldVector< ctype, coorddim > lower)
Constructor from a single point only.
Definition: axisalignedcubegeometry.hh:150
static constexpr int mydimension
Dimension of the cube element.
Definition: axisalignedcubegeometry.hh:56
static constexpr int coorddimension
Dimension of the world space that the cube element is embedded in.
Definition: axisalignedcubegeometry.hh:59
GlobalCoordinate corner(int k) const
Return world coordinates of the k-th corner of the element.
Definition: axisalignedcubegeometry.hh:260
ctype Volume
Type used for volume.
Definition: axisalignedcubegeometry.hh:71
FieldVector< ctype, dim > LocalCoordinate
Type used for a vector of element coordinates.
Definition: axisalignedcubegeometry.hh:65
JacobianTransposed jacobianTransposed(const LocalCoordinate &local) const
Jacobian transposed of the transformation from local to global coordinates.
Definition: axisalignedcubegeometry.hh:196
FieldVector< ctype, coorddim > GlobalCoordinate
Type used for a vector of world coordinates.
Definition: axisalignedcubegeometry.hh:68
LocalCoordinate local(const GlobalCoordinate &global) const
Map a point in global (world) coordinates to element coordinates.
Definition: axisalignedcubegeometry.hh:180
CoordType ctype
Type used for single coordinate coefficients.
Definition: axisalignedcubegeometry.hh:62
std::conditional_t< dim==coorddim, DiagonalMatrix< ctype, dim >, FieldMatrix< ctype, coorddim, dim > > Jacobian
Return type of jacobian.
Definition: axisalignedcubegeometry.hh:100
GeometryType type() const
Type of the cube. Here: a hypercube of the correct dimension.
Definition: axisalignedcubegeometry.hh:155
int corners() const
Return the number of corners of the element.
Definition: axisalignedcubegeometry.hh:254
Jacobian jacobian(const LocalCoordinate &local) const
Jacobian of the transformation from local to global coordinates.
Definition: axisalignedcubegeometry.hh:220
std::conditional< dim==coorddim, DiagonalMatrix< ctype, dim >, FieldMatrix< ctype, coorddim, dim > >::type JacobianInverseTransposed
Return type of jacobianInverseTransposed.
Definition: axisalignedcubegeometry.hh:91
Volume integrationElement(const LocalCoordinate &local) const
Return the integration element, i.e., the determinant term in the integral transformation formula.
Definition: axisalignedcubegeometry.hh:234
GlobalCoordinate center() const
Return center of mass of the element.
Definition: axisalignedcubegeometry.hh:240
bool affine() const
Return if the element is affine. Here: yes.
Definition: axisalignedcubegeometry.hh:302
std::conditional_t< dim==coorddim, DiagonalMatrix< ctype, dim >, FieldMatrix< ctype, dim, coorddim > > JacobianInverse
Return type of jacobianInverse.
Definition: axisalignedcubegeometry.hh:109
GlobalCoordinate global(const LocalCoordinate &local) const
Map a point in local (element) coordinates to world coordinates.
Definition: axisalignedcubegeometry.hh:161
A diagonal matrix of static size.
Definition: diagonalmatrix.hh:53
A dense n x m matrix.
Definition: fmatrix.hh:117
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:126
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.
unspecified-type ReferenceElement
Returns the type of reference element for the argument types T...
Definition: referenceelements.hh:417
constexpr GeometryType cube(unsigned int dim)
Returns a GeometryType representing a hypercube of dimension dim.
Definition: type.hh:473
Mask< V > mask(ADLTag< 0, std::is_same< V, Mask< V > >::value >, const V &v)
implements Simd::mask()
Definition: defaults.hh:153
Dune namespace.
Definition: alignedallocator.hh:13
static const ReferenceElement & cube()
get hypercube reference elements
Definition: referenceelements.hh:210
A unique label for each type of element that can occur in a grid.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 21, 23:30, 2024)