3 #ifndef DUNE_IDENTITYGRIDGEOMETRY_HH
4 #define DUNE_IDENTITYGRIDGEOMETRY_HH
10 #include <dune/common/fmatrix.hh>
11 #include <dune/common/typetraits.hh>
16 template<
int mydim,
int coorddim,
class Gr
idImp>
22 typedef typename GridImp::ctype ctype;
32 typedef typename GridImp::HostGridType::Traits::template Codim<CodimInHostGrid>::Geometry
HostGridGeometryType;
35 typedef typename std::conditional<coorddim==DimensionWorld, HostGridGeometryType, HostGridLocalGeometryType>::type
HostGridGeometry;
67 const FieldVector<ctype, coorddim>
corner (
int i)
const {
74 FieldVector<ctype, coorddim>
global (
const FieldVector<ctype, mydim>& local)
const {
87 FieldVector<ctype, mydim> local (
const FieldVector<ctype, coorddim>&
global)
const {
93 bool checkInside(
const FieldVector<ctype, mydim> &local)
const {
Include standard header files.
Definition: agrid.hh:58
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:130
Default implementation for class Geometry.
Definition: common/geometry.hh:301
Definition: identitygridgeometry.hh:19
bool checkInside(const FieldVector< ctype, mydim > &local) const
Returns true if the point is in the current element.
Definition: identitygridgeometry.hh:93
JacobianInverseTransposed jacobianInverseTransposed(const FieldVector< ctype, mydim > &local) const
The Jacobian matrix of the mapping from the reference element to this element.
Definition: identitygridgeometry.hh:106
JacobianTransposed jacobianTransposed(const FieldVector< ctype, mydim > &local) const
Return the transposed of the Jacobian.
Definition: identitygridgeometry.hh:81
bool affine() const
Definition: identitygridgeometry.hh:56
IdentityGridGeometry(const HostGridGeometry &hostGeometry)
Definition: identitygridgeometry.hh:44
GeometryType type() const
Return the element type identifier.
Definition: identitygridgeometry.hh:51
@ DimensionWorld
Definition: identitygridgeometry.hh:29
FieldVector< ctype, coorddim > global(const FieldVector< ctype, mydim > &local) const
Maps a local coordinate within reference element to global coordinate in element
Definition: identitygridgeometry.hh:74
HostGridGeometryType::JacobianTransposed JacobianTransposed
Definition: identitygridgeometry.hh:39
const FieldVector< ctype, coorddim > corner(int i) const
access to coordinates of corners. Index is the number of the corner
Definition: identitygridgeometry.hh:67
ctype integrationElement(const FieldVector< ctype, mydim > &local) const
Definition: identitygridgeometry.hh:100
HostGridGeometry hostGeometry_
Definition: identitygridgeometry.hh:111
GridImp::HostGridType::Traits::template Codim< CodimInHostGrid >::Geometry HostGridLocalGeometryType
Definition: identitygridgeometry.hh:33
std::conditional< coorddim==DimensionWorld, HostGridGeometryType, HostGridLocalGeometryType >::type HostGridGeometry
Definition: identitygridgeometry.hh:35
GridImp::HostGridType::Traits::template Codim< CodimInHostGrid >::Geometry HostGridGeometryType
Definition: identitygridgeometry.hh:32
@ CodimInHostGrid
Definition: identitygridgeometry.hh:28
HostGridGeometryType::JacobianInverseTransposed JacobianInverseTransposed
type of jacobian transposed
Definition: identitygridgeometry.hh:38
int corners() const
return the number of corners of this element. Corners are numbered 0...n-1
Definition: identitygridgeometry.hh:61
Wrapper and interface classes for element geometries.