This class provides access to geometric and topological properties of a reference element. This includes its type, the number of subentities, the volume, and a method for checking if a point is inside. The embedding of each subentity into the reference element is also provided.
More...
|
| int | size (int c) const |
| | number of subentities of codimension c More...
|
| |
| int | size (int i, int c, int cc) const |
| | number of subentities of codimension cc of subentity (i,c) More...
|
| |
| int | subEntity (int i, int c, int ii, int cc) const |
| | obtain number of ii-th subentity with codim cc of (i,c) More...
|
| |
| const GeometryType & | type (int i, int c) const |
| | obtain the type of subentity (i,c) More...
|
| |
| const GeometryType & | type () const |
| | obtain the type of this reference element More...
|
| |
| const FieldVector< ctype, dim > & | position (int i, int c) const |
| | position of the barycenter of entity (i,c) More...
|
| |
| bool | checkInside (const FieldVector< ctype, dim > &local) const |
| | check if a coordinate is in the reference element More...
|
| |
| template<int codim> |
| bool | checkInside (const FieldVector< ctype, dim-codim > &local, int i) const |
| | check if a local coordinate is in the reference element of the i-th subentity E with codimension c of the current reference element. More...
|
| |
| template<int codim> |
| FieldVector< ctype, dim > | global (const FieldVector< ctype, dim-codim > &local, int i, int c) const |
| | map a local coordinate on subentity (i,codim) into the reference element More...
|
| |
| template<int codim> |
| FieldVector< ctype, dim > | global (const FieldVector< ctype, dim-codim > &local, int i) const |
| | map a local coordinate on subentity (i,codim) into the reference element More...
|
| |
| template<int codim> |
| Codim< codim >::Geometry | geometry (int i) const |
| | obtain the embedding of subentity (i,codim) into the reference element More...
|
| |
| template<int codim> |
| const Codim< codim >::Mapping & | mapping (int i) const |
| | obtain the embedding of subentity (i,codim) into the reference element More...
|
| |
| ctype | volume () const |
| | obtain the volume of the reference element More...
|
| |
| const FieldVector< ctype, dim > & | integrationOuterNormal (int face) const |
| | obtain the integration outer normal of the reference element More...
|
| |
| const FieldVector< ctype, dim > & | volumeOuterNormal (int face) const |
| |
| void | initializeTopology (unsigned int topologyId) |
| | initialize the reference element More...
|
| |
template<class ctype, int dim>
class Dune::ReferenceElement< ctype, dim >
This class provides access to geometric and topological properties of a reference element. This includes its type, the number of subentities, the volume, and a method for checking if a point is inside. The embedding of each subentity into the reference element is also provided.
A singleton of this class for a given geometry type can be accessed through the ReferenceElements class.
- Template Parameters
-
| ctype | field type for coordinates |
| dim | dimension of the reference element |
template<class ctype, int dim>
template<int codim>
| bool Dune::ReferenceElement< ctype, dim >::checkInside |
( |
const FieldVector< ctype, dim-codim > & |
local, |
|
|
int |
i |
|
) |
| const |
|
inline |
check if a local coordinate is in the reference element of the i-th subentity E with codimension c of the current reference element.
Denote by E the i-th subentity of codimension codim of the current reference element. This method return true, if the given local coordinate is within the reference element for the entity E.
- Template Parameters
-
| codim | codimension of subentity E |
- Parameters
-
| [in] | local | coordinates of the point with respect to the reference element of E |
| [in] | i | number of subentity E (0 <= i < size( c )) |
References Dune::GeometryType::id(), and Dune::ReferenceElement< ctype, dim >::type().
template<class ctype, int dim>
template<int codim>
| FieldVector< ctype, dim > Dune::ReferenceElement< ctype, dim >::global |
( |
const FieldVector< ctype, dim-codim > & |
local, |
|
|
int |
i, |
|
|
int |
c |
|
) |
| const |
|
inline |
map a local coordinate on subentity (i,codim) into the reference element
Denote by E the i-th subentity of codimension codim of the current reference element. This method maps a point within the reference element of E into the current reference element.
- Template Parameters
-
| codim | codimension of subentity E |
- Parameters
-
| [in] | local | coordinates of the point with respect to the reference element of E |
| [in] | i | number of subentity E (0 <= i < size( c )) |
| [in] | c | codimension of subentity E |
- Note
- The runtime argument c is redundant and must equal codim.
- Deprecated:
- "Use geometry< codim >( i ).global( local ) instead."
Referenced by Dune::ReferenceElement< ctype, dim >::global().
template<class ctype, int dim>
template<int codim>
| FieldVector< ctype, dim > Dune::ReferenceElement< ctype, dim >::global |
( |
const FieldVector< ctype, dim-codim > & |
local, |
|
|
int |
i |
|
) |
| const |
|
inline |
map a local coordinate on subentity (i,codim) into the reference element
Denote by E the i-th subentity of codimension codim of the current reference element. This method maps a point within the reference element of E into the current reference element.
- Template Parameters
-
| codim | codimension of subentity E |
- Parameters
-
| [in] | local | coordinates of the point with respect to the reference element of E |
| [in] | i | number of subentity E (0 <= i < size( codim )) |
- Deprecated:
- "Use geometry< codim >( i ).global( local ) instead."
References Dune::ReferenceElement< ctype, dim >::global().
template<class ctype, int dim>
position of the barycenter of entity (i,c)
Denote by E the i-th subentity of codimension c of the current reference element. This method returns the coordinates of the center of gravity of E within the current reference element.
- Parameters
-
| [in] | i | number of subentity E (0 <= i < size( c )) |
| [in] | c | codimension of subentity E |
Referenced by Dune::AffineGeometry< ct, mydim, cdim >::center(), Dune::AffineGeometry< ct, mydim, cdim >::corner(), Dune::RefinementImp::Simplex::RefinementIteratorSpecial< dimension, CoordType, 0 >::geometry(), and Dune::MultiLinearGeometry< ct, mydim, cdim, Traits >::local().
template<class ctype, int dim>
number of subentities of codimension cc of subentity (i,c)
Denote by E the i-th subentity of codimension c of the current reference element. This method returns the number of subentities of codimension cc of the current reference element, that are also a subentity of E.
- Parameters
-
| [in] | i | number of subentity E (0 <= i < size( c )) |
| [in] | c | codimension of subentity E |
| [in] | cc | codimension whose size is desired (c <= cc <= dim) |
References Dune::ReferenceElement< ctype, dim >::size().
template<class ctype, int dim>
obtain number of ii-th subentity with codim cc of (i,c)
Denote by E the i-th subentity of codimension c of the current reference element. And denote by S the ii-th subentity of codimension (cc-c) of E. Then, S is a also a subentity of codimension c of the current reference element. This method returns the number of S with respect to the current reference element.
- Parameters
-
| [in] | i | number of subentity E (0 <= i < size( c )) |
| [in] | c | codimension of subentity E |
| [in] | ii | number of subentity S (with respect to E) |
| [in] | cc | codimension of subentity S (c <= cc <= dim) |
References Dune::ReferenceElement< ctype, dim >::size().
template<class ctype, int dim>
obtain the type of subentity (i,c)
Denote by E the i-th subentity of codimension c of the current reference element. This method returns the GeometryType of E.
- Parameters
-
| [in] | i | number of subentity E (0 <= i < size( c )) |
| [in] | c | codimension of subentity E |
References Dune::ReferenceElement< ctype, dim >::size().
Referenced by Dune::RefinementImp::Simplex::RefinementIteratorSpecial< dimension, CoordType, 0 >::geometry(), Dune::MultiLinearGeometry< ct, mydim, cdim, Traits >::topologyId(), and Dune::AffineGeometry< ct, mydim, cdim >::type().