Class holding additionally mpfa data of adaptive compositional models. More...
#include <dumux/porousmediumflow/2p2c/sequential/variableclassadaptive.hh>
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.
TypeTag | The 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... | |
CellDataVector & | cellDataGlobal () |
Return the vector holding all cell data. More... | |
const CellDataVector & | cellDataGlobal () 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, mpfaData > | irregularInterfaceMap_ |
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 |
|
inherited |
|
inline |
Constructs a VariableClass object.
gridView | a DUNE gridview object corresponding to diffusion and transport equation |
|
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.
size | Size of the current (refined and coarsened) grid |
|
inlineinherited |
Return the cell data of a specific cell.
|
inlineinherited |
|
inlineinherited |
Return the vector holding all cell data.
|
inlineinherited |
|
inlineinherited |
Return mapper for elements (for adaptive grids)
|
inlineinherited |
Return mapper for elements (for static grids)
|
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.
irregularIs | The current irregular intersection |
secondHalfEdgeIntersectionIt | Iterator to the intersection connecting the second interaction region |
T1 | Transmissitivity matrix for flux calculations: unique interaction region |
T1_secondHalfEdge | Second transmissitivity matrix for flux calculations for non-unique interaction region |
globalPos3 | The position of the 3rd cell of the first interaction region |
globalIdx3 | The index of the 3rd cell of the first interaction region |
|
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.
irregularIs | The current irregular intersection |
T1 | Transmissitivity matrix for flux calculations: unique interaction region |
globalPos3 | The position of the 3rd cell of the first interaction region |
globalIdx3 | The index of the 3rd cell of the first interaction region |
globalPos4 | The position of the 4th cell of the interaction region |
globalIdx4 | The index of the 4th cell of the interaction region |
subFaceIdx | The index of the subface (up to 4 subfaces possible in 3D) |
|
inlineinherited |
Return gridView.
|
inlineinherited |
Get index of element (codim 0 entity)
Get index of element (codim 0 entity).
element | codim 0 entity |
|
inlineinherited |
Get index of vertex (codim dim entity)
Get index of vertex (codim dim entity).
vertex | codim dim entity |
|
inlineinherited |
Initializes the variable class.
Method initializes the cellData vector. Should be called from problem init()
|
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.
irregularIs | The current irregular intersection |
weight | The weighting factor |
subFaceIdx | The index of the subface (up to 4 subfaces possible in 3D) |
|
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).
problem | The current 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.
problem | The current problem |
|
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.
irregularIs | The current irregular intersection |
T1 | Transmissitivity matrix for flux calculations |
globalPos3 | The position of the 3rd cell of the interaction region |
globalIdx3 | The index of the 3rd cell of the interaction region |
|
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.
irregularIs | The current irregular intersection |
secondHalfEdgeIntersectionIt | Iterator to the intersection connecting the second interaction region |
T1 | Transmissitivity matrix for flux calculations: unique interaction region |
T1_secondHalfEdge | Second transmissitivity matrix for flux calculations for non-unique interaction region |
globalPos3 | The position of the 3rd cell of the first interaction region |
globalIdx3 | The index of the 3rd cell of the first interaction region |
|
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.
irregularIs | The current irregular intersection |
T1 | Transmissitivity matrix for flux calculations |
globalPos3 | The position of the 3rd cell of the interaction region |
globalIdx3 | The index of the 3rd cell of the interaction region |
globalPos4 | The position of the 4th cell of the interaction region |
globalIdx4 | The index of the 4th cell of the interaction region |
subFaceIdx | The index of the subface (up to 4 subfaces possible in 3D) |
|
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.
problem | The current problem |
|
inlineinherited |
Return mapper for vertices (for adaptive grids)
|
inlineinherited |
Return mapper for vertices (for static grids)
|
protected |
The grid.
|
protected |
Storage container for mpfa information.
|
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.