Dune::cpgrid::Geometry< 3, cdim > Class Template Reference

Specialization for 3-dimensional geometries, i.e. cells. More...

#include <Geometry.hpp>

Public Types

enum  { dimension = 3 }
 Dimension of underlying grid.
 
enum  { mydimension = 3 }
 Dimension of domain space of. More...
 
enum  { coorddimension = cdim }
 Dimension of range space of. More...
 
enum  { dimensionworld = 3 }
 World dimension of underlying grid.
 
typedef double ctype
 Coordinate element type.
 
typedef FieldVector< ctype, mydimension > LocalCoordinate
 Domain type of. More...
 
typedef FieldVector< ctype, coorddimension > GlobalCoordinate
 Range type of. More...
 
typedef FieldMatrix< ctype, coorddimension, mydimension > Jacobian
 Type of Jacobian matrix.
 
typedef FieldMatrix< ctype, mydimension, coorddimension > JacobianTransposed
 Type of transposed Jacobian matrix.
 
typedef FieldMatrix< ctype, coorddimension, mydimension > JacobianInverseTransposed
 Type of the inverse of the transposed Jacobian matrix.
 
typedef Dune::GenericGeometry::MatrixHelper< Dune::GenericGeometry::DuneCoordTraits< double > > MatrixHelperType
 

Public Member Functions

 Geometry (const GlobalCoordinate &pos, ctype vol, const GlobalCoordinate *allcorners, const int *corner_indices)
 Construct from centroid, volume (1- and 0-moments) and corners. More...
 
 Geometry (const GlobalCoordinate &pos, ctype vol)
 Construct from centroid and volume (1- and 0-moments). More...
 
 Geometry ()
 Default constructor, giving a non-valid geometry.
 
GlobalCoordinate global (const LocalCoordinate &local_coord) const
 Provide a trilinear mapping. More...
 
LocalCoordinate local (const GlobalCoordinate &y) const
 Mapping from the cell to the reference domain. More...
 
double integrationElement (const LocalCoordinate &local_coord) const
 Equal to {{J^T J}} where J is the Jacobian. More...
 
GeometryType type () const
 Using the cube type for all entities now (cells and vertices), but we use the singular type for intersections. More...
 
int corners () const
 The number of corners of this convex polytope. More...
 
GlobalCoordinate corner (int cor) const
 The 8 corners of the hexahedral base cell.
 
ctype volume () const
 Cell volume.
 
const GlobalCoordinatecenter () const
 Returns the centroid of the geometry.
 
const JacobianTransposed jacobianTransposed (const LocalCoordinate &local_coord) const
 Jacobian transposed. More...
 
const JacobianInverseTransposed jacobianInverseTransposed (const LocalCoordinate &local_coord) const
 Inverse of Jacobian transposed. More...
 
bool affine () const
 The mapping implemented by this geometry is not generally affine.
 

Detailed Description

template<int cdim>
class Dune::cpgrid::Geometry< 3, cdim >

Specialization for 3-dimensional geometries, i.e. cells.

Member Typedef Documentation

◆ GlobalCoordinate

template<int cdim>
typedef FieldVector<ctype, coorddimension> Dune::cpgrid::Geometry< 3, cdim >::GlobalCoordinate

Range type of.

See also
global().

◆ LocalCoordinate

template<int cdim>
typedef FieldVector<ctype, mydimension> Dune::cpgrid::Geometry< 3, cdim >::LocalCoordinate

Domain type of.

See also
global().

Member Enumeration Documentation

◆ anonymous enum

template<int cdim>
anonymous enum

Dimension of domain space of.

See also
global().

◆ anonymous enum

template<int cdim>
anonymous enum

Dimension of range space of.

See also
global().

Constructor & Destructor Documentation

◆ Geometry() [1/2]

template<int cdim>
Dune::cpgrid::Geometry< 3, cdim >::Geometry ( const GlobalCoordinate pos,
ctype  vol,
const GlobalCoordinate allcorners,
const int *  corner_indices 
)
inline

Construct from centroid, volume (1- and 0-moments) and corners.

Parameters
posthe centroid of the entity
volthe volume(area) of the entity
allcornersarray of all corner positions in the grid
corner_indicesarray of 8 indices into allcorners array. The indices must be given in lexicographical order by (kji), i.e. i running fastest.

◆ Geometry() [2/2]

template<int cdim>
Dune::cpgrid::Geometry< 3, cdim >::Geometry ( const GlobalCoordinate pos,
ctype  vol 
)
inline

Construct from centroid and volume (1- and 0-moments).

Note that since corners are not given, the geometry provides no mappings, and some calls (corner(), global() etc.) will fail. This possibly dangerous constructor is available for the benefit of CpGrid::readSintefLegacyFormat().

Parameters
posthe centroid of the entity
volthe volume(area) of the entity

Member Function Documentation

◆ corners()

template<int cdim>
int Dune::cpgrid::Geometry< 3, cdim >::corners ( ) const
inline

The number of corners of this convex polytope.

Returning 8, since we treat all cells as hexahedral.

◆ global()

template<int cdim>
GlobalCoordinate Dune::cpgrid::Geometry< 3, cdim >::global ( const LocalCoordinate local_coord) const
inline

Provide a trilinear mapping.

Note that this does not give a proper space-filling embedding of the cell complex in the general (faulted) case. We should therefore revisit this at some point.

◆ integrationElement()

template<int cdim>
double Dune::cpgrid::Geometry< 3, cdim >::integrationElement ( const LocalCoordinate local_coord) const
inline

Equal to {{J^T J}} where J is the Jacobian.

J_{ij} = (dg_i/du_j) where g is the mapping from the reference domain, and {u_j} are the reference coordinates.

◆ jacobianInverseTransposed()

template<int cdim>
const JacobianInverseTransposed Dune::cpgrid::Geometry< 3, cdim >::jacobianInverseTransposed ( const LocalCoordinate local_coord) const
inline

Inverse of Jacobian transposed.

See also
jacobianTransposed().

◆ jacobianTransposed()

template<int cdim>
const JacobianTransposed Dune::cpgrid::Geometry< 3, cdim >::jacobianTransposed ( const LocalCoordinate local_coord) const
inline

Jacobian transposed.

J^T_{ij} = (dg_j/du_i) where g is the mapping from the reference domain, and {u_i} are the reference coordinates.

◆ local()

template<int cdim>
LocalCoordinate Dune::cpgrid::Geometry< 3, cdim >::local ( const GlobalCoordinate y) const
inline

Mapping from the cell to the reference domain.

May be slow.

◆ type()

template<int cdim>
GeometryType Dune::cpgrid::Geometry< 3, cdim >::type ( ) const
inline

Using the cube type for all entities now (cells and vertices), but we use the singular type for intersections.


The documentation for this class was generated from the following file: