dune-geometry 3.0-git
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
11#include <dune/common/fmatrix.hh>
12#include <dune/common/fvector.hh>
13
14#include <dune/geometry/type.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
52 typedef FieldVector< ctype, mydimension > LocalCoordinate;
53
55 typedef FieldVector< ctype, coorddimension > GlobalCoordinate;
56
58 typedef FieldMatrix< ctype, mydimension, coorddimension > JacobianTransposed;
59
61 typedef FieldMatrix< ctype, coorddimension, mydimension > JacobianInverseTransposed;
62
63 private:
66
68
69 // Helper class to compute a matrix pseudo inverse
71
72 public:
74 AffineGeometry ( const ReferenceElement &refElement, const GlobalCoordinate &origin,
75 const JacobianTransposed &jt )
76 : refElement_(&refElement), origin_(origin), jacobianTransposed_(jt)
77 {
78 integrationElement_ = MatrixHelper::template rightInvA< mydimension, coorddimension >( jacobianTransposed_, jacobianInverseTransposed_ );
79 }
80
83 const JacobianTransposed &jt )
84 : refElement_( &ReferenceElements::general( gt ) ), origin_(origin), jacobianTransposed_( jt )
85 {
86 integrationElement_ = MatrixHelper::template rightInvA< mydimension, coorddimension >( jacobianTransposed_, jacobianInverseTransposed_ );
87 }
88
90 template< class CoordVector >
91 AffineGeometry ( const ReferenceElement &refElement, const CoordVector &coordVector )
92 : refElement_(&refElement), origin_(coordVector[0])
93 {
94 for( int i = 0; i < mydimension; ++i )
95 jacobianTransposed_[ i ] = coordVector[ i+1 ] - origin_;
96 integrationElement_ = MatrixHelper::template rightInvA< mydimension, coorddimension >( jacobianTransposed_, jacobianInverseTransposed_ );
97 }
98
100 template< class CoordVector >
101 AffineGeometry ( Dune::GeometryType gt, const CoordVector &coordVector )
102 : refElement_(&ReferenceElements::general( gt )), origin_(coordVector[0] )
103 {
104 for( int i = 0; i < mydimension; ++i )
105 jacobianTransposed_[ i ] = coordVector[ i+1 ] - origin_;
106 integrationElement_ = MatrixHelper::template rightInvA< mydimension, coorddimension >( jacobianTransposed_, jacobianInverseTransposed_ );
107 }
108
110 bool affine () const { return true; }
111
113 Dune::GeometryType type () const { return refElement_->type(); }
114
116 int corners () const { return refElement_->size( mydimension ); }
117
119 GlobalCoordinate corner ( int i ) const
120 {
121 return global( refElement_->position( i, mydimension ) );
122 }
123
125 GlobalCoordinate center () const { return global( refElement_->position( 0, 0 ) ); }
126
134 {
135 GlobalCoordinate global( origin_ );
136 jacobianTransposed_.umtv( local, global );
137 return global;
138 }
139
154 {
156 jacobianInverseTransposed_.mtv( global - origin_, local );
157 return local;
158 }
159
171 {
172 DUNE_UNUSED_PARAMETER(local);
173 return integrationElement_;
174 }
175
177 ctype volume () const
178 {
179 return integrationElement_ * refElement_->volume();
180 }
181
189 {
190 DUNE_UNUSED_PARAMETER(local);
191 return jacobianTransposed_;
192 }
193
201 {
202 DUNE_UNUSED_PARAMETER(local);
203 return jacobianInverseTransposed_;
204 }
205
206 friend const ReferenceElement &referenceElement ( const AffineGeometry &geometry ) { return *geometry.refElement_; }
207
208 private:
209 const ReferenceElement* refElement_;
210 GlobalCoordinate origin_;
211 JacobianTransposed jacobianTransposed_;
212 JacobianInverseTransposed jacobianInverseTransposed_;
213 ctype integrationElement_;
214 };
215
216} // namespace Dune
217
218#endif // #ifndef DUNE_GEOMETRY_AFFINEGEOMETRY_HH
A unique label for each type of element that can occur in a grid.
Definition affinegeometry.hh:19
This class provides access to geometric and topological properties of a reference element.
Definition referenceelements.hh:53
int size(int c) const
number of subentities of codimension c
Definition referenceelements.hh:80
const FieldVector< ctype, dim > & position(int i, int c) const
position of the barycenter of entity (i,c)
Definition referenceelements.hh:148
ctype volume() const
obtain the volume of the reference element
Definition referenceelements.hh:186
const GeometryType & type(int i, int c) const
obtain the type of subentity (i,c)
Definition referenceelements.hh:130
Class providing access to the singletons of the reference elements.
Definition referenceelements.hh:455
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:91
AffineGeometry(Dune::GeometryType gt, const GlobalCoordinate &origin, const JacobianTransposed &jt)
Create affine geometry from GeometryType, one vertex, and the Jacobian matrix.
Definition affinegeometry.hh:82
friend const ReferenceElement & referenceElement(const AffineGeometry &geometry)
Definition affinegeometry.hh:206
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:113
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:74
ctype volume() const
Obtain the volume of the element.
Definition affinegeometry.hh:177
AffineGeometry(Dune::GeometryType gt, const CoordVector &coordVector)
Create affine geometry from GeometryType and a vector of vertex coordinates.
Definition affinegeometry.hh:101
ctype integrationElement(const LocalCoordinate &local) const
Obtain the integration element.
Definition affinegeometry.hh:170
const JacobianInverseTransposed & jacobianInverseTransposed(const LocalCoordinate &local) const
Obtain the transposed of the Jacobian's inverse.
Definition affinegeometry.hh:200
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:119
int corners() const
Obtain number of corners of the corresponding reference element.
Definition affinegeometry.hh:116
LocalCoordinate local(const GlobalCoordinate &global) const
Evaluate the inverse mapping.
Definition affinegeometry.hh:153
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:133
GlobalCoordinate center() const
Obtain the centroid of the mapping's image.
Definition affinegeometry.hh:125
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
bool affine() const
Always true: this is an affine geometry.
Definition affinegeometry.hh:110
const JacobianTransposed & jacobianTransposed(const LocalCoordinate &local) const
Obtain the transposed of the Jacobian.
Definition affinegeometry.hh:188
Definition matrixhelper.hh:36
Unique label for each type of entities that can occur in DUNE grids.
Definition type.hh:25