Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal > Class Template Reference

Represents the stencil (finite volume geometry) of a single element in the ECFV discretization. More...

#include <ecfvstencil.hh>

Classes

class  EcfvSubControlVolumeFace
 Represents a face of a sub-control volume. More...
 
class  SubControlVolume
 Represents a sub-control volume. More...
 

Public Types

using Entity = Element
 
using Mapper = ElementMapper
 
using LocalGeometry = typename Element::Geometry
 
using SubControlVolumeFace = EcfvSubControlVolumeFace< needFaceIntegrationPos, needFaceNormal >
 
using BoundaryFace = EcfvSubControlVolumeFace< true, needFaceNormal >
 

Public Member Functions

 EcfvStencil (const GridView &gridView, const Mapper &mapper)
 
void updateTopology (const Element &element)
 
void updatePrimaryTopology (const Element &element)
 
void update (const Element &element)
 
void updateCenterGradients ()
 
const Element & element () const
 Return the element to which the stencil refers. More...
 
const GridView & gridView () const
 Return the grid view of the element to which the stencil refers. More...
 
size_t numDof () const
 Returns the number of degrees of freedom which the current element interacts with. More...
 
size_t numPrimaryDof () const
 Returns the number of degrees of freedom which are contained by within the current element. More...
 
unsigned globalSpaceIndex (unsigned dofIdx) const
 Return the global space index given the index of a degree of freedom. More...
 
Dune::PartitionType partitionType (unsigned dofIdx) const
 Return partition type of a given degree of freedom. More...
 
const Element & element (unsigned dofIdx) const
 Return the element given the index of a degree of freedom. More...
 
const Entityentity (unsigned dofIdx) const
 Return the entity given the index of a degree of freedom. More...
 
const SubControlVolumesubControlVolume (unsigned dofIdx) const
 Returns the sub-control volume object belonging to a given degree of freedom. More...
 
size_t numInteriorFaces () const
 Returns the number of interior faces of the stencil. More...
 
const SubControlVolumeFaceinteriorFace (unsigned faceIdx) const
 Returns the face object belonging to a given face index in the interior of the domain. More...
 
size_t numBoundaryFaces () const
 Returns the number of boundary faces of the stencil. More...
 
const BoundaryFaceboundaryFace (unsigned bfIdx) const
 Returns the boundary face object belonging to a given boundary face index. More...
 

Protected Attributes

const GridView & gridView_
 
const ElementMapper & elementMapper_
 
std::vector< Element > elements_
 
std::vector< SubControlVolumesubControlVolumes_
 
std::vector< SubControlVolumeFaceinteriorFaces_
 
std::vector< BoundaryFaceboundaryFaces_
 

Detailed Description

template<class Scalar, class GridView, bool needFaceIntegrationPos = true, bool needFaceNormal = true>
class Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >

Represents the stencil (finite volume geometry) of a single element in the ECFV discretization.

The ECFV discretization is a element centered finite volume approach. This means that each element corresponds to a control volume.

Member Typedef Documentation

◆ BoundaryFace

template<class Scalar , class GridView , bool needFaceIntegrationPos = true, bool needFaceNormal = true>
using Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::BoundaryFace = EcfvSubControlVolumeFace<true, needFaceNormal>

◆ Entity

template<class Scalar , class GridView , bool needFaceIntegrationPos = true, bool needFaceNormal = true>
using Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::Entity = Element

◆ LocalGeometry

template<class Scalar , class GridView , bool needFaceIntegrationPos = true, bool needFaceNormal = true>
using Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::LocalGeometry = typename Element::Geometry

◆ Mapper

template<class Scalar , class GridView , bool needFaceIntegrationPos = true, bool needFaceNormal = true>
using Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::Mapper = ElementMapper

◆ SubControlVolumeFace

template<class Scalar , class GridView , bool needFaceIntegrationPos = true, bool needFaceNormal = true>
using Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::SubControlVolumeFace = EcfvSubControlVolumeFace<needFaceIntegrationPos, needFaceNormal>

Constructor & Destructor Documentation

◆ EcfvStencil()

template<class Scalar , class GridView , bool needFaceIntegrationPos = true, bool needFaceNormal = true>
Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::EcfvStencil ( const GridView &  gridView,
const Mapper mapper 
)
inline

Member Function Documentation

◆ boundaryFace()

template<class Scalar , class GridView , bool needFaceIntegrationPos = true, bool needFaceNormal = true>
const BoundaryFace & Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::boundaryFace ( unsigned  bfIdx) const
inline

Returns the boundary face object belonging to a given boundary face index.

References Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::boundaryFaces_.

◆ element() [1/2]

◆ element() [2/2]

template<class Scalar , class GridView , bool needFaceIntegrationPos = true, bool needFaceNormal = true>
const Element & Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::element ( unsigned  dofIdx) const
inline

Return the element given the index of a degree of freedom.

If no degree of freedom index is passed, the element which was passed to the update() method is returned...

References Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::elements_, and Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::numDof().

◆ entity()

template<class Scalar , class GridView , bool needFaceIntegrationPos = true, bool needFaceNormal = true>
const Entity & Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::entity ( unsigned  dofIdx) const
inline

Return the entity given the index of a degree of freedom.

References Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::element().

◆ globalSpaceIndex()

template<class Scalar , class GridView , bool needFaceIntegrationPos = true, bool needFaceNormal = true>
unsigned Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::globalSpaceIndex ( unsigned  dofIdx) const
inline

◆ gridView()

template<class Scalar , class GridView , bool needFaceIntegrationPos = true, bool needFaceNormal = true>
const GridView & Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::gridView ( ) const
inline

◆ interiorFace()

template<class Scalar , class GridView , bool needFaceIntegrationPos = true, bool needFaceNormal = true>
const SubControlVolumeFace & Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::interiorFace ( unsigned  faceIdx) const
inline

Returns the face object belonging to a given face index in the interior of the domain.

References Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::interiorFaces_.

◆ numBoundaryFaces()

template<class Scalar , class GridView , bool needFaceIntegrationPos = true, bool needFaceNormal = true>
size_t Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::numBoundaryFaces ( ) const
inline

Returns the number of boundary faces of the stencil.

References Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::boundaryFaces_.

◆ numDof()

template<class Scalar , class GridView , bool needFaceIntegrationPos = true, bool needFaceNormal = true>
size_t Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::numDof ( ) const
inline

◆ numInteriorFaces()

template<class Scalar , class GridView , bool needFaceIntegrationPos = true, bool needFaceNormal = true>
size_t Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::numInteriorFaces ( ) const
inline

Returns the number of interior faces of the stencil.

References Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::interiorFaces_.

◆ numPrimaryDof()

template<class Scalar , class GridView , bool needFaceIntegrationPos = true, bool needFaceNormal = true>
size_t Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::numPrimaryDof ( ) const
inline

Returns the number of degrees of freedom which are contained by within the current element.

Primary DOFs are always expected to have a lower index than "secondary" DOFs.

For element centered finite elements, this is only the central DOF.

◆ partitionType()

template<class Scalar , class GridView , bool needFaceIntegrationPos = true, bool needFaceNormal = true>
Dune::PartitionType Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::partitionType ( unsigned  dofIdx) const
inline

Return partition type of a given degree of freedom.

References Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::elements_.

◆ subControlVolume()

template<class Scalar , class GridView , bool needFaceIntegrationPos = true, bool needFaceNormal = true>
const SubControlVolume & Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::subControlVolume ( unsigned  dofIdx) const
inline

Returns the sub-control volume object belonging to a given degree of freedom.

References Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::subControlVolumes_.

◆ update()

template<class Scalar , class GridView , bool needFaceIntegrationPos = true, bool needFaceNormal = true>
void Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::update ( const Element &  element)
inline

◆ updateCenterGradients()

template<class Scalar , class GridView , bool needFaceIntegrationPos = true, bool needFaceNormal = true>
void Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::updateCenterGradients ( )
inline

◆ updatePrimaryTopology()

template<class Scalar , class GridView , bool needFaceIntegrationPos = true, bool needFaceNormal = true>
void Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::updatePrimaryTopology ( const Element &  element)
inline

◆ updateTopology()

Member Data Documentation

◆ boundaryFaces_

template<class Scalar , class GridView , bool needFaceIntegrationPos = true, bool needFaceNormal = true>
std::vector<BoundaryFace> Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::boundaryFaces_
protected

◆ elementMapper_

template<class Scalar , class GridView , bool needFaceIntegrationPos = true, bool needFaceNormal = true>
const ElementMapper& Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::elementMapper_
protected

◆ elements_

◆ gridView_

template<class Scalar , class GridView , bool needFaceIntegrationPos = true, bool needFaceNormal = true>
const GridView& Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::gridView_
protected

◆ interiorFaces_

template<class Scalar , class GridView , bool needFaceIntegrationPos = true, bool needFaceNormal = true>
std::vector<SubControlVolumeFace> Opm::EcfvStencil< Scalar, GridView, needFaceIntegrationPos, needFaceNormal >::interiorFaces_
protected

◆ subControlVolumes_


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