dune-geometry 2.11
Loading...
Searching...
No Matches
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// SPDX-FileCopyrightInfo: Copyright © DUNE Project contributors, see file LICENSE.md in module root
4// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
5#ifndef DUNE_GEOMETRY_AFFINEGEOMETRY_HH
6#define DUNE_GEOMETRY_AFFINEGEOMETRY_HH
7
12
13#include <cmath>
14
15#include <dune/common/fmatrix.hh>
16#include <dune/common/fvector.hh>
17
18#include <dune/geometry/type.hh>
20
21namespace Dune
22{
23
24 // External Forward Declarations
25 // -----------------------------
26
27 namespace Geo
28 {
29
30 template< typename Implementation >
31 class ReferenceElement;
32
33 template< class ctype, int dim >
35
36 template< class ctype, int dim >
37 struct ReferenceElements;
38
39 }
40
41
47 template< class ct, int mydim, int cdim>
49 {
50 public:
51
53 typedef ct ctype;
54
56 static const int mydimension= mydim;
57
59 static const int coorddimension = cdim;
60
62 typedef FieldVector< ctype, mydimension > LocalCoordinate;
63
65 typedef FieldVector< ctype, coorddimension > GlobalCoordinate;
66
68 typedef ctype Volume;
69
71 typedef FieldMatrix< ctype, mydimension, coorddimension > JacobianTransposed;
72
74 typedef FieldMatrix< ctype, coorddimension, mydimension > JacobianInverseTransposed;
75
77 typedef FieldMatrix< ctype, coorddimension, mydimension > Jacobian;
78
80 typedef FieldMatrix< ctype, mydimension, coorddimension > JacobianInverse;
81
82 private:
85
86 typedef Geo::ReferenceElements< ctype, mydimension > ReferenceElements;
87
88 // Helper class to compute a matrix pseudo inverse
89 typedef Impl::FieldMatrixHelper< ct > MatrixHelper;
90
91 public:
98 AffineGeometry () = default;
99
101 AffineGeometry ( const ReferenceElement &refElement, const GlobalCoordinate &origin,
102 const JacobianTransposed &jt )
103 : refElement_(refElement), origin_(origin), jacobianTransposed_(jt)
104 {
105 integrationElement_ = MatrixHelper::rightInvA( jacobianTransposed_, jacobianInverseTransposed_ );
106 }
107
110 const JacobianTransposed &jt )
111 : AffineGeometry(ReferenceElements::general( gt ), origin, jt)
112 { }
113
115 template< class CoordVector >
116 AffineGeometry ( const ReferenceElement &refElement, const CoordVector &coordVector )
117 : refElement_(refElement), origin_(coordVector[0])
118 {
119 for( int i = 0; i < mydimension; ++i )
120 jacobianTransposed_[ i ] = coordVector[ i+1 ] - origin_;
121 integrationElement_ = MatrixHelper::rightInvA( jacobianTransposed_, jacobianInverseTransposed_ );
122 }
123
125 template< class CoordVector >
126 AffineGeometry ( Dune::GeometryType gt, const CoordVector &coordVector )
127 : AffineGeometry(ReferenceElements::general( gt ), coordVector)
128 { }
129
131 bool affine () const { return true; }
132
134 Dune::GeometryType type () const { return refElement_.type(); }
135
137 int corners () const { return refElement_.size( mydimension ); }
138
140 GlobalCoordinate corner ( int i ) const
141 {
142 return global( refElement_.position( i, mydimension ) );
143 }
144
146 GlobalCoordinate center () const { return global( refElement_.position( 0, 0 ) ); }
147
155 {
156 GlobalCoordinate global( origin_ );
157 jacobianTransposed_.umtv( local, global );
158 return global;
159 }
160
175 {
177 jacobianInverseTransposed_.mtv( global - origin_, local );
178 return local;
179 }
180
191 ctype integrationElement ([[maybe_unused]] const LocalCoordinate &local) const
192 {
193 return integrationElement_;
194 }
195
197 Volume volume () const
198 {
199 return integrationElement_ * refElement_.volume();
200 }
201
208 const JacobianTransposed &jacobianTransposed ([[maybe_unused]] const LocalCoordinate &local) const
209 {
210 return jacobianTransposed_;
211 }
212
220 {
221 return jacobianInverseTransposed_;
222 }
223
230 Jacobian jacobian ([[maybe_unused]] const LocalCoordinate &local) const
231 {
232 return jacobianTransposed_.transposed();
233 }
234
241 JacobianInverse jacobianInverse ([[maybe_unused]] const LocalCoordinate &local) const
242 {
243 return jacobianInverseTransposed_.transposed();
244 }
245
246 friend ReferenceElement referenceElement ( const AffineGeometry &geometry )
247 {
248 return geometry.refElement_;
249 }
250
251 private:
252 ReferenceElement refElement_;
253 GlobalCoordinate origin_;
254 JacobianTransposed jacobianTransposed_;
255 JacobianInverseTransposed jacobianInverseTransposed_;
256 ctype integrationElement_;
257 };
258
259} // namespace Dune
260
261#endif // #ifndef DUNE_GEOMETRY_AFFINEGEOMETRY_HH
A unique label for each type of element that can occur in a grid.
Definition affinegeometry.hh:22
decltype(referenceElement(std::declval< T >()...)) ReferenceElement
Definition referenceelements.hh:291
Definition affinegeometry.hh:28
This class provides access to geometric and topological properties of a reference element.
Definition referenceelement.hh:52
Definition affinegeometry.hh:34
Class providing access to the singletons of the reference elements.
Definition referenceelements.hh:128
AffineGeometry(const ReferenceElement &refElement, const CoordVector &coordVector)
Create affine geometry from reference element and a vector of vertex coordinates.
Definition affinegeometry.hh:116
AffineGeometry()=default
Constructs an empty geometry.
AffineGeometry(Dune::GeometryType gt, const GlobalCoordinate &origin, const JacobianTransposed &jt)
Create affine geometry from GeometryType, one vertex, and the Jacobian matrix.
Definition affinegeometry.hh:109
FieldVector< ctype, mydimension > LocalCoordinate
Type for local coordinate vector.
Definition affinegeometry.hh:62
Dune::GeometryType type() const
Obtain the type of the reference element.
Definition affinegeometry.hh:134
static const int mydimension
Dimension of the geometry.
Definition affinegeometry.hh:56
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:101
ctype Volume
Type used for volume.
Definition affinegeometry.hh:68
JacobianInverse jacobianInverse(const LocalCoordinate &local) const
Obtain the Jacobian's inverse.
Definition affinegeometry.hh:241
friend ReferenceElement referenceElement(const AffineGeometry &geometry)
Definition affinegeometry.hh:246
AffineGeometry(Dune::GeometryType gt, const CoordVector &coordVector)
Create affine geometry from GeometryType and a vector of vertex coordinates.
Definition affinegeometry.hh:126
ctype integrationElement(const LocalCoordinate &local) const
Obtain the integration element.
Definition affinegeometry.hh:191
FieldMatrix< ctype, mydimension, coorddimension > JacobianInverse
Type for the inverse Jacobian matrix.
Definition affinegeometry.hh:80
FieldMatrix< ctype, coorddimension, mydimension > Jacobian
Type for the Jacobian matrix.
Definition affinegeometry.hh:77
const JacobianInverseTransposed & jacobianInverseTransposed(const LocalCoordinate &local) const
Obtain the transposed of the Jacobian's inverse.
Definition affinegeometry.hh:219
FieldMatrix< ctype, mydimension, coorddimension > JacobianTransposed
Type for the transposed Jacobian matrix.
Definition affinegeometry.hh:71
GlobalCoordinate corner(int i) const
Obtain coordinates of the i-th corner.
Definition affinegeometry.hh:140
int corners() const
Obtain number of corners of the corresponding reference element.
Definition affinegeometry.hh:137
LocalCoordinate local(const GlobalCoordinate &global) const
Evaluate the inverse mapping.
Definition affinegeometry.hh:174
FieldMatrix< ctype, coorddimension, mydimension > JacobianInverseTransposed
Type for the transposed inverse Jacobian matrix.
Definition affinegeometry.hh:74
static const int coorddimension
Dimension of the world space.
Definition affinegeometry.hh:59
GlobalCoordinate global(const LocalCoordinate &local) const
Evaluate the mapping.
Definition affinegeometry.hh:154
GlobalCoordinate center() const
Obtain the centroid of the mapping's image.
Definition affinegeometry.hh:146
Jacobian jacobian(const LocalCoordinate &local) const
Obtain the Jacobian.
Definition affinegeometry.hh:230
ct ctype
Type used for coordinates.
Definition affinegeometry.hh:53
FieldVector< ctype, coorddimension > GlobalCoordinate
Type for coordinate vector in world space.
Definition affinegeometry.hh:65
bool affine() const
Always true: this is an affine geometry.
Definition affinegeometry.hh:131
const JacobianTransposed & jacobianTransposed(const LocalCoordinate &local) const
Obtain the transposed of the Jacobian.
Definition affinegeometry.hh:208
Volume volume() const
Obtain the volume of the element.
Definition affinegeometry.hh:197
Unique label for each type of entities that can occur in DUNE grids.
Definition type.hh:114