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. More...
 
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. More...
 
typedef double ctype
 Coordinate element type. More...
 
typedef FieldVector< ctype, mydimensionLocalCoordinate
 Domain type of. More...
 
typedef FieldVector< ctype, coorddimensionGlobalCoordinate
 Range type of. More...
 
typedef FieldMatrix< ctype, coorddimension, mydimensionJacobian
 Type of Jacobian matrix. More...
 
typedef FieldMatrix< ctype, coorddimension, mydimensionJacobianInverse
 Type of inverse of Jacobian matrix. More...
 
typedef FieldMatrix< ctype, mydimension, coorddimensionJacobianTransposed
 Type of transposed Jacobian matrix. More...
 
typedef FieldMatrix< ctype, coorddimension, mydimensionJacobianInverseTransposed
 Type of the inverse of the transposed Jacobian matrix. More...
 
typedef Dune::Impl::FieldMatrixHelper< double > MatrixHelperType
 
typedef Dune::FieldVector< double, 3 > PointType
 Refine a single cell considering different widths, lengths, and heights. More...
 

Public Member Functions

 Geometry (const GlobalCoordinate &pos, ctype vol, std::shared_ptr< const EntityVariable< cpgrid::Geometry< 0, 3 >, 3 > > allcorners_ptr, const int *corner_indices)
 Construct from center, volume (1- and 0-moments) and corners. More...
 
 Geometry (const GlobalCoordinate &pos, ctype vol)
 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(). More...
 
 Geometry ()
 Default constructor, giving a non-valid geometry. More...
 
GlobalCoordinate global (const LocalCoordinate &local_coord) const
 
LocalCoordinate local (const GlobalCoordinate &y) const
 
double integrationElement (const LocalCoordinate &local_coord) const
 
GeometryType type () const
 
int corners () const
 
GlobalCoordinate corner (int cor) const
 Get the cor-th of 8 corners of the hexahedral base cell. More...
 
ctype volume () const
 Cell volume. More...
 
void set_volume (ctype volume)
 
const GlobalCoordinatecenter () const
 Returns the centroid of the geometry. More...
 
const JacobianTransposed jacobianTransposed (const LocalCoordinate &local_coord) const
 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. g = g(u) = (g_1(u), g_2(u), g_3(u)), u=(u_1,u_2,u_3) g = map from (local) reference domain to global cell. More...
 
const JacobianInverseTransposed jacobianInverseTransposed (const LocalCoordinate &local_coord) const
 Inverse of Jacobian transposed. More...
 
Jacobian jacobian (const LocalCoordinate &local_coord) const
 The jacobian. More...
 
JacobianInverse jacobianInverse (const LocalCoordinate &local_coord) const
 The inverse of the jacobian. More...
 
bool affine () const
 The mapping implemented by this geometry is not generally affine. More...
 
void refineCellifiedPatch (const std::array< int, 3 > &cells_per_dim, DefaultGeometryPolicy &all_geom, std::vector< std::array< int, 8 > > &refined_cell_to_point, cpgrid::OrientedEntityTable< 0, 1 > &refined_cell_to_face, Opm::SparseTable< int > &refined_face_to_point, cpgrid::OrientedEntityTable< 1, 0 > &refined_face_to_cell, cpgrid::EntityVariable< enum face_tag, 1 > &refined_face_tags, cpgrid::SignedEntityVariable< PointType, 1 > &refined_face_normals, const std::array< int, 3 > &patch_dim, const std::vector< double > &widthsX, const std::vector< double > &lengthsY, const std::vector< double > &heightsZ) const
 

Detailed Description

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

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

Member Typedef Documentation

◆ ctype

template<int cdim>
typedef double Dune::cpgrid::Geometry< 3, cdim >::ctype

Coordinate element type.

◆ GlobalCoordinate

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

Range type of.

See also
global().

◆ Jacobian

template<int cdim>
typedef FieldMatrix< ctype, coorddimension, mydimension > Dune::cpgrid::Geometry< 3, cdim >::Jacobian

Type of Jacobian matrix.

◆ JacobianInverse

template<int cdim>
typedef FieldMatrix< ctype, coorddimension, mydimension > Dune::cpgrid::Geometry< 3, cdim >::JacobianInverse

Type of inverse of Jacobian matrix.

◆ JacobianInverseTransposed

template<int cdim>
typedef FieldMatrix< ctype, coorddimension, mydimension > Dune::cpgrid::Geometry< 3, cdim >::JacobianInverseTransposed

Type of the inverse of the transposed Jacobian matrix.

◆ JacobianTransposed

template<int cdim>
typedef FieldMatrix< ctype, mydimension, coorddimension > Dune::cpgrid::Geometry< 3, cdim >::JacobianTransposed

Type of transposed Jacobian matrix.

◆ LocalCoordinate

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

Domain type of.

See also
global().

◆ MatrixHelperType

template<int cdim>
typedef Dune::Impl::FieldMatrixHelper< double > Dune::cpgrid::Geometry< 3, cdim >::MatrixHelperType

◆ PointType

template<int cdim>
typedef Dune::FieldVector<double,3> Dune::cpgrid::Geometry< 3, cdim >::PointType

Refine a single cell considering different widths, lengths, and heights.

For each cell to be created, storage must be passed for its corners and the indices. That storage must be externally managed, since the newly created geometry structures only store pointers and do not free them on destruction.

Parameters
cells_per_dimThe number of sub-cells in each direction,
all_geomGeometry Policy for the refined geometries. Those will be added there.
refined_cell_to_pointMap from cell to its 8 corners.
refined_cell_to_faceMap from cell to its oriented faces, used to build face_to_cell_.
refined_face_to_pointMap from face to its points.
refined_face_to_cellMap from face to its neighboring cells.
refined_face_tagsFace tags (I_FACE, J_FACE, K_FACE).
refined_face_normalsFace normal(s) (only one per face is computed).
patch_dimAmount of cells to be refined, in each direction.
dx,dy,dzVectors of widths (x-dir), lengths (y-dir), and heights (z-dir)

Member Enumeration Documentation

◆ anonymous enum

template<int cdim>
anonymous enum

Dimension of underlying grid.

Enumerator
dimension 

◆ anonymous enum

template<int cdim>
anonymous enum

Dimension of domain space of.

See also
global().
Enumerator
mydimension 

◆ anonymous enum

template<int cdim>
anonymous enum

Dimension of range space of.

See also
global().
Enumerator
coorddimension 

◆ anonymous enum

template<int cdim>
anonymous enum

World dimension of underlying grid.

Enumerator
dimensionworld 

Constructor & Destructor Documentation

◆ Geometry() [1/3]

template<int cdim>
Dune::cpgrid::Geometry< 3, cdim >::Geometry ( const GlobalCoordinate pos,
ctype  vol,
std::shared_ptr< const EntityVariable< cpgrid::Geometry< 0, 3 >, 3 > >  allcorners_ptr,
const int *  corner_indices 
)
inline

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

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

◆ Geometry() [2/3]

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

◆ Geometry() [3/3]

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

Default constructor, giving a non-valid geometry.

Member Function Documentation

◆ affine()

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

The mapping implemented by this geometry is not generally affine.

◆ center()

template<int cdim>
const GlobalCoordinate & Dune::cpgrid::Geometry< 3, cdim >::center ( ) const
inline

Returns the centroid of the geometry.

◆ corner()

template<int cdim>
GlobalCoordinate Dune::cpgrid::Geometry< 3, cdim >::corner ( int  cor) const
inline

Get the cor-th of 8 corners of the hexahedral base cell.

◆ 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. Map g from (local) reference domain to (global) cell

◆ integrationElement()

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

Equal to \sqrt{\det{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.

◆ jacobian()

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

The jacobian.

◆ jacobianInverse()

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

The inverse of the jacobian.

◆ 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. g = g(u) = (g_1(u), g_2(u), g_3(u)), u=(u_1,u_2,u_3) g = map from (local) reference domain to global cell.

◆ 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.

◆ refineCellifiedPatch()

template<int cdim>
void Dune::cpgrid::Geometry< 3, cdim >::refineCellifiedPatch ( const std::array< int, 3 > &  cells_per_dim,
DefaultGeometryPolicy all_geom,
std::vector< std::array< int, 8 > > &  refined_cell_to_point,
cpgrid::OrientedEntityTable< 0, 1 > &  refined_cell_to_face,
Opm::SparseTable< int > &  refined_face_to_point,
cpgrid::OrientedEntityTable< 1, 0 > &  refined_face_to_cell,
cpgrid::EntityVariable< enum face_tag, 1 > &  refined_face_tags,
cpgrid::SignedEntityVariable< PointType, 1 > &  refined_face_normals,
const std::array< int, 3 > &  patch_dim,
const std::vector< double > &  widthsX,
const std::vector< double > &  lengthsY,
const std::vector< double > &  heightsZ 
) const
inline

— REFINED CORNERS —

— END REFINED CORNERS —

— REFINED FACES —

— END REFINED FACES —

— REFINED CELLS —

— END REFINED CELLS —

References Dune::cpgrid::OrientedEntityTable< codim_from, codim_to >::appendRow(), Opm::SparseTable< T >::appendRow(), Dune::area(), Dune::cpgrid::DefaultGeometryPolicy::geomVector(), J_FACE, Dune::simplex_volume(), and Dune::volume().

◆ set_volume()

template<int cdim>
void Dune::cpgrid::Geometry< 3, cdim >::set_volume ( ctype  volume)
inline

References Dune::volume().

◆ 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.

◆ volume()

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

Cell volume.


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