|
opm-upscaling
|
Specialization for 2D quadrilaterals. More...
#include <boundarygrid.hh>
Public Types | |
| enum | { dimension = 2 } |
| The dimension of the grid. | |
| enum | { mydimension = 2 } |
| Dimension of the domain space. | |
| enum | { coorddimension = cdim } |
| Dimension of the range space. | |
| enum | { dimensionworld = 2 } |
| World dimension of underlying grid. | |
| typedef double | ctype |
| Coordinate element type. | |
| typedef Dune::FieldVector< ctype, mydimension > | LocalCoordinate |
| Domain type. | |
| typedef Dune::FieldVector< ctype, coorddimension > | GlobalCoordinate |
| Range type. | |
| typedef Dune::FieldMatrix< ctype, coorddimension, mydimension > | Jacobian |
| Type of Jacobian matrix. | |
| typedef Dune::FieldMatrix< ctype, mydimension, coorddimension > | JacobianTransposed |
| Type of transposed Jacobian matrix. | |
Public Member Functions | |
| HexGeometry (const BoundaryGrid::Quad &q, const GridImp &gv, int dir) | |
| Construct integration element extracted from a 3D grid. | |
| HexGeometry (const BoundaryGrid::Quad &q) | |
| Construct integration element. | |
| Dune::GeometryType | type () const |
| Returns entity type (a 2D cube). | |
| int | corners () const |
| Returns number of corners. | |
| ctype | volume () const |
| Returns volume (area) of quadrilateral. | |
| GlobalCoordinate | center () const |
| Returns center of quadrilateral. | |
| GlobalCoordinate | corner (int cor) const |
| Returns coordinates to requested corner. | |
| GlobalCoordinate | global (const LocalCoordinate &local) const |
| Map from local coordinates to global coordinates. | |
| LocalCoordinate | local (const GlobalCoordinate &y) const |
| Map from global coordinates to local coordinates. | |
| const Dune::FieldMatrix< ctype, mydimension, coorddimension > | jacobianTransposed (const LocalCoordinate &local) const |
| Return the transposed jacobian. | |
| const Dune::FieldMatrix< ctype, coorddimension, mydimension > | jacobianInverseTransposed (const LocalCoordinate &local) const |
| Returns the inverse, transposed Jacobian. | |
| ctype | integrationElement (const LocalCoordinate &local) const |
| Returns the integration element (|J'*J|)^(1/2). | |
Specialization for 2D quadrilaterals.
|
inline |
Construct integration element extracted from a 3D grid.
| [in] | q | Quad describing element |
| [in] | gv | Underlying 3D grid quads are extracted from |
| [in] | dir | The direction of the normal vector on the face |
|
inlineexplicit |
Construct integration element.
| [in] | q | Quad describing element |
|
inline |
Returns coordinates to requested corner.
| [in] | cor | The requested corner (0..3) |
|
inline |
Map from local coordinates to global coordinates.
| [in] | local | The local coordinates |
|
inline |
Returns the integration element (|J'*J|)^(1/2).
| [in] | local | The local coordinates |
|
inline |
Returns the inverse, transposed Jacobian.
| [in] | local | The local coordinates |
|
inline |
Return the transposed jacobian.
| [in] | local | The local coordinates |
|
inline |
Map from global coordinates to local coordinates.
| [in] | y | The global coordinates |