Class including the information of a 3d interaction volume of an adaptive MPFA L-method that does not change with time. More...
#include <dumux/porousmediumflow/sequential/cellcentered/mpfa/linteractionvolume3dadaptive.hh>

Class including the information of a 3d interaction volume of an adaptive MPFA L-method that does not change with time.
Includes information needed to calculate the transmissibility matrices of an L-interaction-volume.
Public Types | |
| enum | HangingNodeTypes { noHangingNode = -1 , twoSmallCells = 0 , fourSmallCellsFace = 1 , fourSmallCellsEdge = 2 , fourSmallCellsDiag = 3 , sixSmallCells = 4 } |
| The different hanging node interaction volume types (see dissertation M. Wolff, http://elib.uni-stuttgart.de/opus/volltexte/2013/8661/) More... | |
| enum | FaceTypes { inside = 1 , boundary = 0 , outside = -1 } |
| enum | { subVolumeTotalNum = IndexTranslator::subVolumeTotalNum , fluxFacesTotalNum = IndexTranslator::fluxFacesTotalNum , fluxEdgesTotalNum = IndexTranslator::fluxEdgesTotalNum } |
Public Member Functions | |
| FvMpfaL3dInteractionVolumeAdaptive (const Grid &grid) | |
| Constructs a FvMpfaL3dInteractionVolumeAdaptive object. More... | |
| void | reset () |
| Reset the interaction volume (deletes stored data) More... | |
| void | setSubVolumeElement (const Element &element, int subVolumeIdx) |
| Store a dune element as a sub volume element. More... | |
| void | setHangingNodeType (int hNType) |
| Store the type of hanging-node-interaction volume. More... | |
| bool | sameLevel () |
| Check if elements in the interaction volume are of the same grid level. More... | |
| bool | hasLevel (int level) |
| Check if an element of a certain grid level is stored. More... | |
| bool | isHangingNodeVolume () |
| Check whether the interaction volume is a hanging-node volume. More... | |
| int | getHangingNodeType () |
| The type of the interaction volume as type of FvMpfaL3dInteractionVolumeAdaptive<TypeTag>::HangingNodeTypes. More... | |
| void | printInteractionVolumeInfo () |
| Print the stored interaction volume data. More... | |
| void | setCenterPosition (const DimVector ¢erVertexPos) |
| Store position of the central vertex. More... | |
| void | setFacePosition (const DimVector &pos, int fIdx) |
| Store the position of a flux face. More... | |
| void | setEdgePosition (const DimVector &pos, int edgeIdx) |
| Store the center of the edges. More... | |
| void | setFaceArea (Scalar faceArea, int fIdx) |
| Store the flux face areas. More... | |
| void | setNormal (DimVector &normal, int subVolumeIdx, int subVolumeFaceIdx) |
| Store the normals. More... | |
| void | setBoundary (BoundaryTypes &boundaryTypes, int subVolumeFaceIdx) |
| Store the types of boundary conditions. More... | |
| void | setOutsideFace (int subVolumeFaceIdx) |
| Define a flux face to be outside the model domain (for boundary vertices) More... | |
| void | setDirichletCondition (PrimaryVariables &condition, int subVolumeFaceIdx) |
| Store Dirichlet boundary conditions. More... | |
| void | setNeumannCondition (PrimaryVariables &condition, int subVolumeFaceIdx) |
| Store Neumann boundary conditions. More... | |
| void | setIndexOnElement (int indexInInside, int subVolumeIdx, int subVolumeFaceIdx) |
| Store the local dune face index dependent on the local interaction volume indices. More... | |
| DimVector & | getCenterPosition () |
| The position of the central vertex. More... | |
| int | getFaceIndexFromSubVolume (int subVolumeIdx, int subVolumeFaceIdx) |
| The local interaction volume face index dependent on the local sub volume indices. More... | |
| int | getEdgeIndexFromSubVolumeFace (int subVolumeIdx, int subVolumeFaceIdx, int subVolumeEdgeIdx) |
| The local interaction volume edge index dependent on the local sub volume indices. More... | |
| int | getElementNumber () |
| Number of stored dune elements. More... | |
| int | getIndexOnElement (int subVolumeIdx, int subVolumeFaceIdx) |
| The local dune face index dependent on the local interaction volume indices. More... | |
| Element | getSubVolumeElement (int subVolumeIdx) |
| The dune element of the interaction volume sub-volume. More... | |
| bool | hasSubVolumeElement (int subVolumeIdx) |
| Check if a dune element is stored for an interaction volume sub-volume. More... | |
| BoundaryTypes & | getBoundaryType (int subVolumeFaceIdx) |
| The boundary types. More... | |
| bool | isOutsideFace (int subVolumeFaceIdx) |
| Check if an interaction volume face is outside the model domain (for boundary vertices) More... | |
| bool | isInsideFace (int subVolumeFaceIdx) |
| Check if an interaction volume face is inside the model domain (for boundary vertices) More... | |
| bool | isBoundaryFace (int subVolumeFaceIdx) |
| Check if an interaction volume face is a model domain boundary (for boundary vertices) More... | |
| bool | isInnerVolume () |
| Check if an interaction volume is not located around a boundary vertex. More... | |
| bool | isBoundaryInteractionVolume () |
| Check if an interaction volume is located around a boundary vertex. More... | |
| PrimaryVariables & | getDirichletValues (int subVolumeFaceIdx) |
| The Dirichlet boundary values. More... | |
| PrimaryVariables & | getNeumannValues (int subVolumeFaceIdx) |
| The Neumann boundary values. More... | |
| DimVector & | getNormal (int subVolumeIdx, int subVolumeFaceIdx) |
| The face normal. More... | |
| DimVector & | getFacePosition (int subVolumeIdx, int subVolumeFaceIdx) |
| The position of the face center. More... | |
| DimVector & | getFacePosition (int fIdx) |
| The position of the face center. More... | |
| DimVector & | getEdgePosition (int subVolumeIdx, int subVolumeFaceIdx, int subVolumeEdgeIdx) |
| The position of the edge center. More... | |
| DimVector & | getEdgePosition (int edgeIdx) |
| The position of the edge center. More... | |
| Scalar & | getFaceArea (int subVolumeIdx, int subVolumeFaceIdx) |
| The interaction volume flux face area. More... | |
| Scalar & | getFaceArea (int fIdx) |
| The interaction volume flux face area. More... | |
|
inherited |
|
inherited |
| enum Dumux::FvMpfaL3dInteractionVolumeAdaptive::HangingNodeTypes |
The different hanging node interaction volume types (see dissertation M. Wolff, http://elib.uni-stuttgart.de/opus/volltexte/2013/8661/)
|
inline |
Constructs a FvMpfaL3dInteractionVolumeAdaptive object.
|
inlineinherited |
The boundary types.
| subVolumeFaceIdx | The local face index in the interaction volume |
|
inlineinherited |
The position of the central vertex.
|
inlineinherited |
The Dirichlet boundary values.
| subVolumeFaceIdx | The local face index in the interaction volume |
|
inlineinherited |
The local interaction volume edge index dependent on the local sub volume indices.
|
inlineinherited |
The position of the edge center.
| edgeIdx | The local edge index in the interaction volume |
|
inlineinherited |
The position of the edge center.
| subVolumeIdx | The local element index in the interaction volume |
| subVolumeFaceIdx | The local face index in the interaction volume element |
| subVolumeEdgeIdx | The local edge index in the interaction volume element |
|
inlineinherited |
Number of stored dune elements.
|
inlineinherited |
The interaction volume flux face area.
| fIdx | The local face index in the interaction volume |
|
inlineinherited |
The interaction volume flux face area.
| subVolumeIdx | The local element index in the interaction volume |
| subVolumeFaceIdx | The local face index in the interaction volume element |
|
inlineinherited |
The local interaction volume face index dependent on the local sub volume indices.
|
inlineinherited |
The position of the face center.
| fIdx | The local face index in the interaction volume |
|
inlineinherited |
The position of the face center.
| subVolumeIdx | The local element index in the interaction volume |
| subVolumeFaceIdx | The local face index in the interaction volume element |
|
inline |
The type of the interaction volume as type of FvMpfaL3dInteractionVolumeAdaptive<TypeTag>::HangingNodeTypes.
|
inlineinherited |
The local dune face index dependent on the local interaction volume indices.
| subVolumeIdx | The local element index in the interaction volume |
| subVolumeFaceIdx | The local face index in the interaction volume element |
|
inlineinherited |
The Neumann boundary values.
| subVolumeFaceIdx | The local face index in the interaction volume |
|
inlineinherited |
The face normal.
| subVolumeIdx | The local element index in the interaction volume |
| subVolumeFaceIdx | The local face index in the interaction volume element |
|
inlineinherited |
The dune element of the interaction volume sub-volume.
| subVolumeIdx | The local element index in the interaction volume |
|
inline |
Check if an element of a certain grid level is stored.
| level | the dune grid level |
true if an element of a certain grid level is stored
|
inlineinherited |
Check if a dune element is stored for an interaction volume sub-volume.
| subVolumeIdx | The local element index in the interaction volume |
true if an element pointer is stored at position subVolumeIdx.
|
inlineinherited |
Check if an interaction volume face is a model domain boundary (for boundary vertices)
| subVolumeFaceIdx | The local face index in the interaction volume |
true if the flux face is a boundary face.
|
inlineinherited |
Check if an interaction volume is located around a boundary vertex.
true if an interaction volume located around a boundary vertex
|
inline |
Check whether the interaction volume is a hanging-node volume.
true if the interaction volume is a hanging-node volume
|
inlineinherited |
Check if an interaction volume is not located around a boundary vertex.
true if an interaction volume is not located around a boundary vertex
|
inlineinherited |
Check if an interaction volume face is inside the model domain (for boundary vertices)
| subVolumeFaceIdx | The local face index in the interaction volume |
true if the flux face is inside the model domain.
|
inlineinherited |
Check if an interaction volume face is outside the model domain (for boundary vertices)
| subVolumeFaceIdx | The local face index in the interaction volume |
true if the flux face is outside the model domain.
|
inline |
Print the stored interaction volume data.
|
inline |
Reset the interaction volume (deletes stored data)
|
inline |
Check if elements in the interaction volume are of the same grid level.
true if all elements are on the same grid level
|
inlineinherited |
Store the types of boundary conditions.
| boundaryTypes | BoundaryTypes object |
| subVolumeFaceIdx | The local face index in the interaction volume |
|
inlineinherited |
Store position of the central vertex.
|
inlineinherited |
Store Dirichlet boundary conditions.
| condition | Vector of primary variables |
| subVolumeFaceIdx | The local face index in the interaction volume |
|
inlineinherited |
Store the center of the edges.
| pos | Position of face center |
| edgeIdx | The interaction volume edge index |
|
inlineinherited |
Store the flux face areas.
| faceArea | The flux face area |
| fIdx | The interaction volume face index |
|
inlineinherited |
Store the position of a flux face.
| pos | Position of face center |
| fIdx | The interaction volume face index |
|
inline |
Store the type of hanging-node-interaction volume.
| hNType | the type of hanging-node-interaction volume of type FvMpfaL3dInteractionVolumeAdaptive<TypeTag>::HangingNodeTypes |
|
inlineinherited |
Store the local dune face index dependent on the local interaction volume indices.
| indexInInside | Face index of the Dune reference element |
| subVolumeIdx | The local element index in the interaction volume |
| subVolumeFaceIdx | The local face index in the interaction volume element |
|
inlineinherited |
Store Neumann boundary conditions.
| condition | Vector phase fluxes |
| subVolumeFaceIdx | The local face index in the interaction volume |
|
inlineinherited |
Store the normals.
| normal | The normal vector |
| subVolumeIdx | The local element index in the interaction volume |
| subVolumeFaceIdx | The local face index in the interaction volume element |
|
inlineinherited |
Define a flux face to be outside the model domain (for boundary vertices)
| subVolumeFaceIdx | The local face index in the interaction volume |
|
inline |
Store a dune element as a sub volume element.
| element | The element |
| subVolumeIdx | The local element index in the interaction volume |