Dune Core Modules (2.3.1)

affinegeometry.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#ifndef DUNE_GEOMETRY_AFFINEGEOMETRY_HH
4#define DUNE_GEOMETRY_AFFINEGEOMETRY_HH
5
13
14#include <dune/geometry/type.hh>
15#include <dune/geometry/genericgeometry/geometrytraits.hh>
16#include <dune/geometry/genericgeometry/matrixhelper.hh>
17
18namespace Dune
19{
20
21 // External Forward Declarations
22 // -----------------------------
23
24 template< class ctype, int dim >
25 class ReferenceElement;
26
27 template< class ctype, int dim >
28 struct ReferenceElements;
29
30
31
37 template< class ct, int mydim, int cdim>
39 {
40 public:
41
43 typedef ct ctype;
44
46 static const int mydimension= mydim;
47
49 static const int coorddimension = cdim;
50
53
56
59
62
65
66 private:
69
71
72 // Helper class to compute a matrix pseudo inverse
73 typedef GenericGeometry::MatrixHelper< GenericGeometry::DuneCoordTraits< ct > > MatrixHelper;
74
75 public:
77 AffineGeometry ( const ReferenceElement &refElement, const GlobalCoordinate &origin,
78 const JacobianTransposed &jt )
79 : refElement_(&refElement), origin_(origin), jacobianTransposed_(jt)
80 {
81 integrationElement_ = MatrixHelper::template rightInvA< mydimension, coorddimension >( jacobianTransposed_, jacobianInverseTransposed_ );
82 }
83
86 const JacobianTransposed &jt )
87 : refElement_( &ReferenceElements::general( gt ) ), origin_(origin), jacobianTransposed_( jt )
88 {
89 integrationElement_ = MatrixHelper::template rightInvA< mydimension, coorddimension >( jacobianTransposed_, jacobianInverseTransposed_ );
90 }
91
93 template< class CoordVector >
94 AffineGeometry ( const ReferenceElement &refElement, const CoordVector &coordVector )
95 : refElement_(&refElement), origin_(coordVector[0])
96 {
97 for( int i = 0; i < mydimension; ++i )
98 jacobianTransposed_[ i ] = coordVector[ i+1 ] - origin_;
99 integrationElement_ = MatrixHelper::template rightInvA< mydimension, coorddimension >( jacobianTransposed_, jacobianInverseTransposed_ );
100 }
101
103 template< class CoordVector >
104 AffineGeometry ( Dune::GeometryType gt, const CoordVector &coordVector )
105 : refElement_(&ReferenceElements::general( gt )), origin_(coordVector[0] )
106 {
107 for( int i = 0; i < mydimension; ++i )
108 jacobianTransposed_[ i ] = coordVector[ i+1 ] - origin_;
109 integrationElement_ = MatrixHelper::template rightInvA< mydimension, coorddimension >( jacobianTransposed_, jacobianInverseTransposed_ );
110 }
111
113 bool affine () const { return true; }
114
116 Dune::GeometryType type () const { return refElement_->type(); }
117
119 int corners () const { return refElement_->size( mydimension ); }
120
122 GlobalCoordinate corner ( int i ) const
123 {
124 return global( refElement_->position( i, mydimension ) );
125 }
126
128 GlobalCoordinate center () const { return global( refElement_->position( 0, 0 ) ); }
129
137 {
138 GlobalCoordinate global( origin_ );
139 jacobianTransposed_.umtv( local, global );
140 return global;
141 }
142
155 {
156 LocalCoordinate local;
157 jacobianInverseTransposed_.mtv( global - origin_, local );
158 return local;
159 }
160
172 {
174 return integrationElement_;
175 }
176
178 ctype volume () const
179 {
180 return integrationElement_ * refElement_->volume();
181 }
182
190 {
192 return jacobianTransposed_;
193 }
194
202 {
204 return jacobianInverseTransposed_;
205 }
206
207 private:
208 const ReferenceElement* refElement_;
209 GlobalCoordinate origin_;
210 JacobianTransposed jacobianTransposed_;
211 JacobianInverseTransposed jacobianInverseTransposed_;
212 ctype integrationElement_;
213
214 };
215
216} // namespace Dune
217
218#endif // #ifndef DUNE_GEOMETRY_AFFINEGEOMETRY_HH
Implementation of the Geometry interface for affine geometries.
Definition: affinegeometry.hh:39
AffineGeometry(const ReferenceElement &refElement, const CoordVector &coordVector)
Create affine geometry from reference element and a vector of vertex coordinates.
Definition: affinegeometry.hh:94
AffineGeometry(Dune::GeometryType gt, const GlobalCoordinate &origin, const JacobianTransposed &jt)
Create affine geometry from GeometryType, one vertex, and the Jacobian matrix.
Definition: affinegeometry.hh:85
FieldVector< ctype, mydimension > LocalCoordinate
Type for local coordinate vector.
Definition: affinegeometry.hh:52
Dune::GeometryType type() const
Obtain the type of the reference element.
Definition: affinegeometry.hh:116
static const int mydimension
Dimension of the geometry.
Definition: affinegeometry.hh:46
AffineGeometry(const ReferenceElement &refElement, const GlobalCoordinate &origin, const JacobianTransposed &jt)
Create affine geometry from reference element, one vertex, and the Jacobian matrix.
Definition: affinegeometry.hh:77
ctype volume() const
Obtain the volume of the element.
Definition: affinegeometry.hh:178
AffineGeometry(Dune::GeometryType gt, const CoordVector &coordVector)
Create affine geometry from GeometryType and a vector of vertex coordinates.
Definition: affinegeometry.hh:104
ctype integrationElement(const LocalCoordinate &local) const
Obtain the integration element.
Definition: affinegeometry.hh:171
const JacobianInverseTransposed & jacobianInverseTransposed(const LocalCoordinate &local) const
Obtain the transposed of the Jacobian's inverse.
Definition: affinegeometry.hh:201
FieldMatrix< ctype, mydimension, coorddimension > JacobianTransposed
Type for the transposed Jacobian matrix.
Definition: affinegeometry.hh:58
GlobalCoordinate corner(int i) const
Obtain coordinates of the i-th corner.
Definition: affinegeometry.hh:122
int corners() const
Obtain number of corners of the corresponding reference element.
Definition: affinegeometry.hh:119
FieldMatrix< ctype, coorddimension, mydimension > JacobianInverseTransposed
Type for the transposed inverse Jacobian matrix.
Definition: affinegeometry.hh:61
static const int coorddimension
Dimension of the world space.
Definition: affinegeometry.hh:49
GlobalCoordinate global(const LocalCoordinate &local) const
Evaluate the mapping.
Definition: affinegeometry.hh:136
GlobalCoordinate center() const
Obtain the centroid of the mapping's image.
Definition: affinegeometry.hh:128
ct ctype
Type used for coordinates.
Definition: affinegeometry.hh:43
FieldVector< ctype, coorddimension > GlobalCoordinate
Type for coordinate vector in world space.
Definition: affinegeometry.hh:55
JacobianInverseTransposed Jacobian
Definition: affinegeometry.hh:64
bool affine() const
Always true: this is an affine geometry.
Definition: affinegeometry.hh:113
const JacobianTransposed & jacobianTransposed(const LocalCoordinate &local) const
Obtain the transposed of the Jacobian.
Definition: affinegeometry.hh:189
void mtv(const X &x, Y &y) const
y = A^T x
Definition: densematrix.hh:413
void umtv(const X &x, Y &y) const
y += A^T x
Definition: densematrix.hh:450
vector space out of a tensor product of fields.
Definition: fvector.hh:92
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:25
This class provides access to geometric and topological properties of a reference element....
Definition: referenceelements.hh:58
int size(int c) const
number of subentities of codimension c
Definition: referenceelements.hh:86
const FieldVector< ctype, dim > & position(int i, int c) const
position of the barycenter of entity (i,c)
Definition: referenceelements.hh:154
ctype volume() const
obtain the volume of the reference element
Definition: referenceelements.hh:280
const GeometryType & type(int i, int c) const
obtain the type of subentity (i,c)
Definition: referenceelements.hh:136
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.
bool gt(const T &first, const T &second, typename EpsilonType< T >::Type epsilon)
test if first greater than second
Definition: float_cmp.cc:132
Dune namespace.
Definition: alignment.hh:14
Class providing access to the singletons of the reference elements. Special methods are available for...
Definition: referenceelements.hh:563
A unique label for each type of element that can occur in a grid.
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentional unused function parameters with.
Definition: unused.hh:18
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 12, 23:30, 2024)