5#ifndef DUNE_IDENTITYGRIDGEOMETRY_HH
6#define DUNE_IDENTITYGRIDGEOMETRY_HH
12#include <dune/common/fmatrix.hh>
13#include <dune/common/typetraits.hh>
18 template<
int mydim,
int coorddim,
class Gr
idImp>
24 typedef typename GridImp::ctype ctype;
34 typedef typename GridImp::HostGridType::Traits::template Codim<CodimInHostGrid>::Geometry
HostGridGeometryType;
37 typedef typename std::conditional<coorddim==DimensionWorld, HostGridGeometryType, HostGridLocalGeometryType>::type
HostGridGeometry;
69 const FieldVector<ctype, coorddim>
corner (
int i)
const {
76 FieldVector<ctype, coorddim>
global (
const FieldVector<ctype, mydim>& local)
const {
89 FieldVector<ctype, mydim> local (
const FieldVector<ctype, coorddim>&
global)
const {
95 bool checkInside(
const FieldVector<ctype, mydim> &local)
const {
Include standard header files.
Definition: agrid.hh:60
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:132
Default implementation for class Geometry.
Definition: common/geometry.hh:406
Definition: identitygridgeometry.hh:21
bool checkInside(const FieldVector< ctype, mydim > &local) const
Returns true if the point is in the current element.
Definition: identitygridgeometry.hh:95
JacobianInverseTransposed jacobianInverseTransposed(const FieldVector< ctype, mydim > &local) const
The Jacobian matrix of the mapping from the reference element to this element.
Definition: identitygridgeometry.hh:108
JacobianTransposed jacobianTransposed(const FieldVector< ctype, mydim > &local) const
Return the transposed of the Jacobian.
Definition: identitygridgeometry.hh:83
std::conditional< coorddim==DimensionWorld, HostGridGeometryType, HostGridLocalGeometryType >::type HostGridGeometry
Definition: identitygridgeometry.hh:37
bool affine() const
Definition: identitygridgeometry.hh:58
static constexpr int CodimInHostGrid
Definition: identitygridgeometry.hh:30
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:76
IdentityGridGeometry(const HostGridGeometry &hostGeometry)
Definition: identitygridgeometry.hh:46
GeometryType type() const
Return the element type identifier.
Definition: identitygridgeometry.hh:53
const FieldVector< ctype, coorddim > corner(int i) const
access to coordinates of corners. Index is the number of the corner
Definition: identitygridgeometry.hh:69
HostGridGeometryType::JacobianTransposed JacobianTransposed
Definition: identitygridgeometry.hh:41
ctype integrationElement(const FieldVector< ctype, mydim > &local) const
Definition: identitygridgeometry.hh:102
HostGridGeometry hostGeometry_
Definition: identitygridgeometry.hh:113
GridImp::HostGridType::Traits::template Codim< CodimInHostGrid >::Geometry HostGridLocalGeometryType
Definition: identitygridgeometry.hh:35
GridImp::HostGridType::Traits::template Codim< CodimInHostGrid >::Geometry HostGridGeometryType
Definition: identitygridgeometry.hh:34
HostGridGeometryType::JacobianInverseTransposed JacobianInverseTransposed
type of jacobian transposed
Definition: identitygridgeometry.hh:40
static constexpr int DimensionWorld
Definition: identitygridgeometry.hh:31
int corners() const
return the number of corners of this element. Corners are numbered 0...n-1
Definition: identitygridgeometry.hh:63
Wrapper and interface classes for element geometries.