3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
Classes | Public Types | Public Member Functions | Protected Attributes | Static Protected Attributes | List of all members
Dumux::VariableClass2P2CAdaptive< TypeTag > Class Template Reference

Class holding additionally mpfa data of adaptive compositional models. More...

#include <dumux/porousmediumflow/2p2c/sequential/variableclassadaptive.hh>

Inheritance diagram for Dumux::VariableClass2P2CAdaptive< TypeTag >:
Inheritance graph

Description

template<class TypeTag>
class Dumux::VariableClass2P2CAdaptive< TypeTag >

Class holding additionally mpfa data of adaptive compositional models.

This class provides the possibility to store and load the transmissibilities (and associated infos) of the mpfa method per irregular face. This class provides the storage container and access methods for both 2D and 3D implementation. Note that according to the number of half-edges (sub-faces) regarded, one ore more transmissibility can be stored and loaded.

Template Parameters
TypeTagThe Type Tag

Classes

struct  mpfaData
 Storage object for data related to the MPFA method. More...
 

Public Types

using CellDataVector = typename std::vector< CellData >
 

Public Member Functions

 VariableClass2P2CAdaptive (const GridView &gridView)
 Constructs a VariableClass object. More...
 
void adaptVariableSize (int size)
 Resizes sequential variable vectors. More...
 
void reconstructPrimVars (Problem &problem)
 Reconstruct missing primary variables (where elements are created/deleted) More...
 
void storeMpfaData (const typename GridView::Intersection &irregularIs, const TransmissivityMatrix &T1, const GlobalPosition &globalPos3, const int &globalIdx3)
 Stores Mpfa Data of one interaction region on an intersection. More...
 
void storeMpfaData (const typename GridView::Intersection &irregularIs, IntersectionIterator &secondHalfEdgeIntersectionIt, const TransmissivityMatrix &T1, const TransmissivityMatrix &T1_secondHalfEdge, const GlobalPosition &globalPos3, const int &globalIdx3)
 Stores Mpfa Data on an intersection for both half-edges. More...
 
void storeMpfaData3D (const typename GridView::Intersection &irregularIs, const TransmissivityMatrix &T1, const GlobalPosition &globalPos3, const int &globalIdx3, const GlobalPosition &globalPos4, const int &globalIdx4, int subFaceIdx=0)
 Stores 3D Mpfa Data on an intersection. More...
 
void performTransmissitivityWeighting (const typename GridView::Intersection &irregularIs, Scalar weight, int subFaceIdx=-1)
 Weigths the transmissivity coefficient by the flux area (3D) More...
 
int getMpfaData (const Intersection &irregularIs, IntersectionIterator &secondHalfEdgeIntersectionIt, TransmissivityMatrix &T1, TransmissivityMatrix &T1_secondHalfEdge, GlobalPosition &globalPos3, int &globalIdx3)
 Provides access to stored Mpfa Data on an intersection for both half-edges. More...
 
int getMpfaData3D (const Intersection &irregularIs, TransmissivityMatrix &T1, GlobalPosition &globalPos3, int &globalIdx3, GlobalPosition &globalPos4, int &globalIdx4, int subFaceIdx=0)
 Provides access to stored 3D Mpfa Data on an intersection for up to 4 subfaces. More...
 
void storePrimVars (const Problem &problem)
 
void reconstructPrimVars (const Problem &problem)
 
void initialize ()
 Initializes the variable class. More...
 
CellDataVectorcellDataGlobal ()
 Return the vector holding all cell data. More...
 
const CellDataVectorcellDataGlobal () const
 
CellData & cellData (const int idx)
 Return the cell data of a specific cell. More...
 
const CellData & cellData (const int idx) const
 
int index (const Element &element) const
 Get index of element (codim 0 entity) More...
 
int index (const Vertex &vertex) const
 Get index of vertex (codim dim entity) More...
 
const GridView & gridView () const
 Return gridView. More...
 
ElementMapper & elementMapper ()
 Return mapper for elements (for adaptive grids) More...
 
const ElementMapper & elementMapper () const
 Return mapper for elements (for static grids) More...
 
VertexMapper & vertexMapper ()
 Return mapper for vertices (for adaptive grids) More...
 
const VertexMapper & vertexMapper () const
 Return mapper for vertices (for static grids) More...
 

Protected Attributes

std::map< IdType, mpfaDatairregularInterfaceMap_
 Storage container for mpfa information. More...
 
const Grid & grid_
 The grid. More...
 

Static Protected Attributes

static const int storageRequirement = Dune::StaticPower<dim-1, dim>::power
 

Member Typedef Documentation

◆ CellDataVector

template<class TypeTag >
using Dumux::VariableClass< TypeTag >::CellDataVector = typename std::vector<CellData>
inherited

Constructor & Destructor Documentation

◆ VariableClass2P2CAdaptive()

template<class TypeTag >
Dumux::VariableClass2P2CAdaptive< TypeTag >::VariableClass2P2CAdaptive ( const GridView &  gridView)
inline

Constructs a VariableClass object.

Parameters
gridViewa DUNE gridview object corresponding to diffusion and transport equation

Member Function Documentation

◆ adaptVariableSize()

template<class TypeTag >
void Dumux::VariableClass2P2CAdaptive< TypeTag >::adaptVariableSize ( int  size)
inline

Resizes sequential variable vectors.

Method that change the size of the vectors for h-adaptive simulations, and clears recently stored transmissibility matrices for the newly adapted grid.

Parameters
sizeSize of the current (refined and coarsened) grid

◆ cellData() [1/2]

template<class TypeTag >
CellData & Dumux::VariableClass< TypeTag >::cellData ( const int  idx)
inlineinherited

Return the cell data of a specific cell.

◆ cellData() [2/2]

template<class TypeTag >
const CellData & Dumux::VariableClass< TypeTag >::cellData ( const int  idx) const
inlineinherited

◆ cellDataGlobal() [1/2]

template<class TypeTag >
CellDataVector & Dumux::VariableClass< TypeTag >::cellDataGlobal ( )
inlineinherited

Return the vector holding all cell data.

◆ cellDataGlobal() [2/2]

template<class TypeTag >
const CellDataVector & Dumux::VariableClass< TypeTag >::cellDataGlobal ( ) const
inlineinherited

◆ elementMapper() [1/2]

template<class TypeTag >
ElementMapper & Dumux::VariableClass< TypeTag >::elementMapper ( )
inlineinherited

Return mapper for elements (for adaptive grids)

◆ elementMapper() [2/2]

template<class TypeTag >
const ElementMapper & Dumux::VariableClass< TypeTag >::elementMapper ( ) const
inlineinherited

Return mapper for elements (for static grids)

◆ getMpfaData()

template<class TypeTag >
int Dumux::VariableClass2P2CAdaptive< TypeTag >::getMpfaData ( const Intersection &  irregularIs,
IntersectionIterator &  secondHalfEdgeIntersectionIt,
TransmissivityMatrix &  T1,
TransmissivityMatrix &  T1_secondHalfEdge,
GlobalPosition &  globalPos3,
int &  globalIdx3 
)
inline

Provides access to stored Mpfa Data on an intersection for both half-edges.

The method gets information of both interaction regions (Transmissitivity as well as details of the 3rd cells of both regions) from a storage container. The key to each element is the index of the intersection, seen from the smaller cell (only this is unique). If we arrive from the "wrong" (i.e. non-unique) direction, we invert fluxes.

Parameters
irregularIsThe current irregular intersection
secondHalfEdgeIntersectionItIterator to the intersection connecting the second interaction region
T1Transmissitivity matrix for flux calculations: unique interaction region
T1_secondHalfEdgeSecond transmissitivity matrix for flux calculations for non-unique interaction region
globalPos3The position of the 3rd cell of the first interaction region
globalIdx3The index of the 3rd cell of the first interaction region

◆ getMpfaData3D()

template<class TypeTag >
int Dumux::VariableClass2P2CAdaptive< TypeTag >::getMpfaData3D ( const Intersection &  irregularIs,
TransmissivityMatrix &  T1,
GlobalPosition &  globalPos3,
int &  globalIdx3,
GlobalPosition &  globalPos4,
int &  globalIdx4,
int  subFaceIdx = 0 
)
inline

Provides access to stored 3D Mpfa Data on an intersection for up to 4 subfaces.

The method gets information of up to 4 interaction regions (Transmissitivity as well as details of the 3rd & 4th cells of the regions) from a storage container. The key to each element is the index of the intersection, seen from the smaller cell (only this is unique). If we arrive from the "wrong" (i.e. non-unique) direction, we invert fluxes.

Parameters
irregularIsThe current irregular intersection
T1Transmissitivity matrix for flux calculations: unique interaction region
globalPos3The position of the 3rd cell of the first interaction region
globalIdx3The index of the 3rd cell of the first interaction region
globalPos4The position of the 4th cell of the interaction region
globalIdx4The index of the 4th cell of the interaction region
subFaceIdxThe index of the subface (up to 4 subfaces possible in 3D)

◆ gridView()

template<class TypeTag >
const GridView & Dumux::VariableClass< TypeTag >::gridView ( ) const
inlineinherited

Return gridView.

◆ index() [1/2]

template<class TypeTag >
int Dumux::VariableClass< TypeTag >::index ( const Element &  element) const
inlineinherited

Get index of element (codim 0 entity)

Get index of element (codim 0 entity).

Parameters
elementcodim 0 entity
Returns
element index

◆ index() [2/2]

template<class TypeTag >
int Dumux::VariableClass< TypeTag >::index ( const Vertex &  vertex) const
inlineinherited

Get index of vertex (codim dim entity)

Get index of vertex (codim dim entity).

Parameters
vertexcodim dim entity
Returns
vertex index

◆ initialize()

template<class TypeTag >
void Dumux::VariableClass< TypeTag >::initialize ( )
inlineinherited

Initializes the variable class.

Method initializes the cellData vector. Should be called from problem init()

◆ performTransmissitivityWeighting()

template<class TypeTag >
void Dumux::VariableClass2P2CAdaptive< TypeTag >::performTransmissitivityWeighting ( const typename GridView::Intersection &  irregularIs,
Scalar  weight,
int  subFaceIdx = -1 
)
inline

Weigths the transmissivity coefficient by the flux area (3D)

Each cell face could contain up to 4 interaction regions in 3D. If fewer interaction regions are considered (e.g. at a boundary), the flux goes through a higher area than was regarded in the calculation of the interaction region. Therefore, the transmissibility coefficients have to be weighted (scaled). If no specific face Index is given, all present Transmissibility matrices are weighted homogeneously by the given weight.

Parameters
irregularIsThe current irregular intersection
weightThe weighting factor
subFaceIdxThe index of the subface (up to 4 subfaces possible in 3D)

◆ reconstructPrimVars() [1/2]

template<class TypeTag >
void Dumux::VariableClassAdaptive< TypeTag >::reconstructPrimVars ( const Problem &  problem)
inlineinherited

Reconstruct missing primary variables (where elements are created/deleted)

To reconstruct the solution in father elements, problem properties might need to be accessed. Starting from the lowest level, the old solution is mapped on the new grid: Where coarsened, new cells get information from old father element. Where refined, a new solution is reconstructed from the old father cell, and then a new son is created. That is then stored into the general data structure (CellData).

Parameters
problemThe current problem

◆ reconstructPrimVars() [2/2]

template<class TypeTag >
void Dumux::VariableClass2P2CAdaptive< TypeTag >::reconstructPrimVars ( Problem &  problem)
inline

Reconstruct missing primary variables (where elements are created/deleted)

To reconstruct the solution in father elements, problem properties might need to be accessed.

Parameters
problemThe current problem

◆ storeMpfaData() [1/2]

template<class TypeTag >
void Dumux::VariableClass2P2CAdaptive< TypeTag >::storeMpfaData ( const typename GridView::Intersection &  irregularIs,
const TransmissivityMatrix &  T1,
const GlobalPosition &  globalPos3,
const int &  globalIdx3 
)
inline

Stores Mpfa Data of one interaction region on an intersection.

The method stores information of ONE interaction region (Transmissitivity as well as details of the 3rd cell in the region) into a storage container. The key to each element is the index of the intersection, seen from the smaller cell (only this is unique). If we arrive from the "wrong" (i.e. non-unique) direction, we invert fluxes.

Parameters
irregularIsThe current irregular intersection
T1Transmissitivity matrix for flux calculations
globalPos3The position of the 3rd cell of the interaction region
globalIdx3The index of the 3rd cell of the interaction region

◆ storeMpfaData() [2/2]

template<class TypeTag >
void Dumux::VariableClass2P2CAdaptive< TypeTag >::storeMpfaData ( const typename GridView::Intersection &  irregularIs,
IntersectionIterator &  secondHalfEdgeIntersectionIt,
const TransmissivityMatrix &  T1,
const TransmissivityMatrix &  T1_secondHalfEdge,
const GlobalPosition &  globalPos3,
const int &  globalIdx3 
)
inline

Stores Mpfa Data on an intersection for both half-edges.

The method stores information of both interaction regions (Transmissitivity as well as details of the 3rd cells of both regions) into a storage container. The key to each element is the index of the intersection, seen from the smaller cell (only this is unique). If we arrive from the "wrong" (i.e. non-unique) direction, we invert fluxes.

Parameters
irregularIsThe current irregular intersection
secondHalfEdgeIntersectionItIterator to the intersection connecting the second interaction region
T1Transmissitivity matrix for flux calculations: unique interaction region
T1_secondHalfEdgeSecond transmissitivity matrix for flux calculations for non-unique interaction region
globalPos3The position of the 3rd cell of the first interaction region
globalIdx3The index of the 3rd cell of the first interaction region

◆ storeMpfaData3D()

template<class TypeTag >
void Dumux::VariableClass2P2CAdaptive< TypeTag >::storeMpfaData3D ( const typename GridView::Intersection &  irregularIs,
const TransmissivityMatrix &  T1,
const GlobalPosition &  globalPos3,
const int &  globalIdx3,
const GlobalPosition &  globalPos4,
const int &  globalIdx4,
int  subFaceIdx = 0 
)
inline

Stores 3D Mpfa Data on an intersection.

The method stores information of the interaction region (Transmissitivity as well as details of the 3rd & 4th cell in the region) into a storage container. The key to each element is the index of the intersection, seen from the smaller cell (only this index is unique). If we arrive from the "wrong" (i.e. non-unique) direction, we invert fluxes.

Parameters
irregularIsThe current irregular intersection
T1Transmissitivity matrix for flux calculations
globalPos3The position of the 3rd cell of the interaction region
globalIdx3The index of the 3rd cell of the interaction region
globalPos4The position of the 4th cell of the interaction region
globalIdx4The index of the 4th cell of the interaction region
subFaceIdxThe index of the subface (up to 4 subfaces possible in 3D)

◆ storePrimVars()

template<class TypeTag >
void Dumux::VariableClassAdaptive< TypeTag >::storePrimVars ( const Problem &  problem)
inlineinherited

Store primary variables

To reconstruct the solution in father elements, problem properties might need to be accessed. From upper level on downwards, the old solution is stored into an container object, before the grid is adapted. Father elements hold averaged information from the son cells for the case of the sons being coarsened.

Parameters
problemThe current problem

◆ vertexMapper() [1/2]

template<class TypeTag >
VertexMapper & Dumux::VariableClass< TypeTag >::vertexMapper ( )
inlineinherited

Return mapper for vertices (for adaptive grids)

◆ vertexMapper() [2/2]

template<class TypeTag >
const VertexMapper & Dumux::VariableClass< TypeTag >::vertexMapper ( ) const
inlineinherited

Return mapper for vertices (for static grids)

Member Data Documentation

◆ grid_

template<class TypeTag >
const Grid& Dumux::VariableClass2P2CAdaptive< TypeTag >::grid_
protected

The grid.

◆ irregularInterfaceMap_

template<class TypeTag >
std::map<IdType, mpfaData> Dumux::VariableClass2P2CAdaptive< TypeTag >::irregularInterfaceMap_
protected

Storage container for mpfa information.

◆ storageRequirement

template<class TypeTag >
const int Dumux::VariableClass2P2CAdaptive< TypeTag >::storageRequirement = Dune::StaticPower<dim-1, dim>::power
staticprotected

in the 2D case, we need to store 1 additional cell. In 3D, we store dim-1=2 cells for one interaction volume, and 8 cells if four interaction volumes are regarded.


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