Base class for all staggered finite-volume problems. More...
#include <dumux/common/staggeredfvproblem.hh>
Base class for all staggered finite-volume problems.
Public Types | |
using | SpatialParams = GetPropType< TypeTag, Properties::SpatialParams > |
Public Member Functions | |
StaggeredFVProblem (std::shared_ptr< const GridGeometry > gridGeometry, const std::string ¶mGroup="") | |
Constructor. More... | |
bool | isDirichletCell (const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolume &scv, int pvIdx) const |
Returns whether a fixed Dirichlet value shall be used at a given cell. More... | |
template<class ElementVolumeVariables , class ElementFaceVariables , class Entity > | |
NumEqVector | source (const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elementFaceVars, const Entity &e) const |
Evaluate the source term for all phases within a given sub-control-volume (-face). More... | |
NumEqVector | neumann (const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolumeFace &scvf) const |
Evaluate the boundary conditions for a neumann boundary segment. More... | |
template<class Entity > | |
PrimaryVariables | initial (const Entity &entity) const |
Evaluate the initial value for an element (for cell-centered primary variables) or face (for velocities) More... | |
template<class SolutionVector > | |
void | applyInitialSolution (SolutionVector &sol) const |
Applies the initial solution for all degrees of freedom of the grid. More... | |
template<class SolutionVector > | |
void | applyInitialCellCenterSolution (SolutionVector &sol, const SubControlVolume &scv, const PrimaryVariables &initSol) const |
Applies the initial cell center solution. More... | |
template<class SolutionVector > | |
void | applyInitialFaceSolution (SolutionVector &sol, const SubControlVolumeFace &scvf, const PrimaryVariables &initSol) const |
Applies the initial face solution. More... | |
const SpatialParams & | spatialParams () const |
Return a reference to the underlying spatial parameters. More... | |
SpatialParams & | spatialParams () |
Return a reference to the underlying spatial parameters. More... | |
const std::string & | name () const |
The problem name. More... | |
void | setName (const std::string &newName) |
Set the problem name. More... | |
Boundary conditions and sources defining the problem | |
auto | boundaryTypes (const Element &element, const SubControlVolume &scv) const |
Specifies which kind of boundary condition should be used for which equation on a given boundary segment. More... | |
auto | boundaryTypes (const Element &element, const SubControlVolumeFace &scvf) const |
Specifies which kind of boundary condition should be used for which equation on a given boundary segment. More... | |
BoundaryTypes | boundaryTypesAtPos (const GlobalPosition &globalPos) const |
Specifies which kind of boundary condition should be used for which equation on a given boundary segment. More... | |
PrimaryVariables | dirichlet (const Element &element, const SubControlVolumeFace &scvf) const |
Evaluate the boundary conditions for a dirichlet control volume face. More... | |
PrimaryVariables | dirichlet (const Element &element, const SubControlVolume &scv) const |
Evaluate the boundary conditions for a dirichlet control volume. More... | |
PrimaryVariables | dirichletAtPos (const GlobalPosition &globalPos) const |
Evaluate the boundary conditions for a dirichlet control volume. More... | |
template<class ElementVolumeVariables , class ElementFluxVariablesCache > | |
NumEqVector | neumann (const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const |
Evaluate the boundary conditions for a neumann boundary segment. More... | |
NumEqVector | neumannAtPos (const GlobalPosition &globalPos) const |
Evaluate the boundary conditions for a neumann boundary segment. More... | |
template<class ElementVolumeVariables > | |
NumEqVector | source (const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv) const |
Evaluate the source term for all phases within a given sub-control-volume. More... | |
NumEqVector | sourceAtPos (const GlobalPosition &globalPos) const |
Evaluate the source term for all phases within a given sub-control-volume. More... | |
void | addPointSources (std::vector< PointSource > &pointSources) const |
Applies a vector of point sources. The point sources are possibly solution dependent. More... | |
template<class ElementVolumeVariables > | |
void | pointSource (PointSource &source, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv) const |
Evaluate the point sources (added by addPointSources) for all phases within a given sub-control-volume. More... | |
void | pointSourceAtPos (PointSource &pointSource, const GlobalPosition &globalPos) const |
Evaluate the point sources (added by addPointSources) for all phases within a given sub-control-volume. More... | |
template<class MatrixBlock , class VolumeVariables > | |
void | addSourceDerivatives (MatrixBlock &block, const Element &element, const FVElementGeometry &fvGeometry, const VolumeVariables &volVars, const SubControlVolume &scv) const |
Add source term derivative to the Jacobian. More... | |
template<class ElementVolumeVariables > | |
NumEqVector | scvPointSources (const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv) const |
Adds contribution of point sources for a specific sub control volume to the values. Caution: Only overload this method in the implementation if you know what you are doing. More... | |
void | computePointSourceMap () |
Compute the point source map, i.e. which scvs have point source contributions. More... | |
const PointSourceMap & | pointSourceMap () const |
Get the point source map. It stores the point sources per scv. More... | |
PrimaryVariables | initialAtPos (const GlobalPosition &globalPos) const |
Evaluate the initial value for a control volume. More... | |
const GridGeometry & | gridGeometry () const |
The finite volume grid geometry. More... | |
const std::string & | paramGroup () const |
The parameter group in which to retrieve runtime parameters. More... | |
static constexpr bool | enableInternalDirichletConstraints () |
If internal Dirichlet constraints are enabled Enables / disables internal (non-boundary) Dirichlet constraints. If this is overloaded to return true, the assembler calls problem.hasInternalDirichletConstraint(element, scv). This means you have to implement the following member function. More... | |
Implementation & | asImp_ () |
Returns the implementation of the problem (i.e. static polymorphism) More... | |
const Implementation & | asImp_ () const |
Returns the implementation of the problem (i.e. static polymorphism) More... | |
|
inherited |
|
inline |
Constructor.
gridGeometry | The finite volume grid geometry |
paramGroup | The parameter group in which to look for runtime parameters first (default is "") |
|
inlineinherited |
Applies a vector of point sources. The point sources are possibly solution dependent.
pointSources | A vector of PointSource s that contain source values for all phases and space positions. |
For this method, the values method of the point source has to return the absolute rate values in units \( [ \textnormal{unit of conserved quantity} / s ] \). Positive values mean that the conserved quantity is created, negative ones mean that it vanishes. E.g. for the mass balance that would be a mass rate in \( [ kg / s ] \).
|
inlineinherited |
Add source term derivative to the Jacobian.
|
inline |
Applies the initial cell center solution.
|
inline |
Applies the initial face solution.
|
inline |
Applies the initial solution for all degrees of freedom of the grid.
|
inlineprotectedinherited |
Returns the implementation of the problem (i.e. static polymorphism)
|
inlineprotectedinherited |
Returns the implementation of the problem (i.e. static polymorphism)
|
inlineinherited |
Specifies which kind of boundary condition should be used for which equation on a given boundary segment.
element | The finite element |
scv | The sub control volume |
|
inlineinherited |
Specifies which kind of boundary condition should be used for which equation on a given boundary segment.
element | The finite element |
scvf | The sub control volume face |
|
inlineinherited |
Specifies which kind of boundary condition should be used for which equation on a given boundary segment.
globalPos | The position of the finite volume in global coordinates |
As a default, i.e. if the user's problem does not overload any boundaryTypes method set Dirichlet boundary conditions everywhere for all primary variables
|
inlineinherited |
Compute the point source map, i.e. which scvs have point source contributions.
|
inlineinherited |
Evaluate the boundary conditions for a dirichlet control volume.
element | The finite element |
scv | the sub control volume |
|
inlineinherited |
Evaluate the boundary conditions for a dirichlet control volume face.
element | The finite element |
scvf | the sub control volume face |
|
inlineinherited |
Evaluate the boundary conditions for a dirichlet control volume.
globalPos | The position of the center of the finite volume for which the dirichlet condition ought to be set in global coordinates |
|
inlinestaticconstexprinherited |
If internal Dirichlet constraints are enabled Enables / disables internal (non-boundary) Dirichlet constraints. If this is overloaded to return true, the assembler calls problem.hasInternalDirichletConstraint(element, scv). This means you have to implement the following member function.
bool hasInternalDirichletConstraint(const Element& element, const SubControlVolume& scv) const;
which returns an indexable container of booleans defining for each equation if the corresponding dof associated with the element/scv pair is constraint. If true is returned for a dof, the assembler calls problem.internalDirichlet(element, scv). This means you have to additionally implement the following member function
PrimaryVariables internalDirichlet(const Element& element, const SubControlVolume& scv) const;
which returns the enforced Dirichlet values the dof associated with the element/scv pair.
|
inlineinherited |
The finite volume grid geometry.
|
inline |
Evaluate the initial value for an element (for cell-centered primary variables) or face (for velocities)
entity | The dof entity (element or vertex) |
|
inlineinherited |
Evaluate the initial value for a control volume.
globalPos | The global position |
|
inline |
Returns whether a fixed Dirichlet value shall be used at a given cell.
element | The finite element |
fvGeometry | The finite-volume geometry |
scv | The sub control volume |
pvIdx | The primary variable index |
|
inlineinherited |
The problem name.
This is used as a prefix for files generated by the simulation. It could be either overwritten by the problem files, or simply declared over the setName() function in the application file.
|
inline |
Evaluate the boundary conditions for a neumann boundary segment.
This is the method for the case where the Neumann condition is potentially solution dependent
element | The finite element |
fvGeometry | The finite-volume geometry |
elemVolVars | All volume variables for the element |
elemFaceVars | All face variables for the element |
scvf | The sub control volume face |
Negative values mean influx. E.g. for the mass balance that would the mass flux in \( [ kg / (m^2 \cdot s)] \).
|
inlineinherited |
Evaluate the boundary conditions for a neumann boundary segment.
This is the method for the case where the Neumann condition is potentially solution dependent
element | The finite element |
fvGeometry | The finite-volume geometry |
elemVolVars | All volume variables for the element |
elemFluxVarsCache | Flux variables caches for all faces in stencil |
scvf | The sub control volume face |
Negative values mean influx. E.g. for the mass balance that would be the mass flux in \( [ kg / (m^2 \cdot s)] \).
|
inlineinherited |
Evaluate the boundary conditions for a neumann boundary segment.
globalPos | The position of the boundary face's integration point in global coordinates |
Negative values mean influx. E.g. for the mass balance that would be the mass flux in \( [ kg / (m^2 \cdot s)] \).
As a default, i.e. if the user's problem does not overload any neumann method return no-flow Neumann boundary conditions at all Neumann boundaries
|
inlineinherited |
The parameter group in which to retrieve runtime parameters.
|
inlineinherited |
Evaluate the point sources (added by addPointSources) for all phases within a given sub-control-volume.
This is the method for the case where the point source is solution dependent
source | A single point source |
element | The finite element |
fvGeometry | The finite-volume geometry |
elemVolVars | All volume variables for the element |
scv | The sub control volume |
For this method, the values() method of the point sources returns the absolute conserved quantity rate generated or annihilate in units \( [ \textnormal{unit of conserved quantity} / s ] \). Positive values mean that the conserved quantity is created, negative ones mean that it vanishes. E.g. for the mass balance that would be a mass rate in \( [ kg / s ] \).
|
inlineinherited |
Evaluate the point sources (added by addPointSources) for all phases within a given sub-control-volume.
This is the method for the case where the point source is space dependent
pointSource | A single point source |
globalPos | The point source position in global coordinates |
For this method, the values() method of the point sources returns the absolute conserved quantity rate generated or annihilate in units \( [ \textnormal{unit of conserved quantity} / s ] \). Positive values mean that the conserved quantity is created, negative ones mean that it vanishes. E.g. for the mass balance that would be a mass rate in \( [ kg / s ] \).
|
inlineinherited |
Get the point source map. It stores the point sources per scv.
|
inlineinherited |
Adds contribution of point sources for a specific sub control volume to the values. Caution: Only overload this method in the implementation if you know what you are doing.
|
inlineinherited |
Set the problem name.
This static method sets the simulation name, which should be called before the application problem is declared! If not, the default name "sim" will be used.
newName | The problem's name |
|
inline |
Evaluate the source term for all phases within a given sub-control-volume (-face).
This is the method for the case where the source term is potentially solution dependent and requires some quantities that are specific to the fully-implicit method.
element | The finite element |
fvGeometry | The finite-volume geometry |
elemVolVars | All volume variables for the element |
elementFaceVars | All face variables for the element |
e | The geometrical entity on which the source shall be applied (scv or scvf) |
For this method, the return parameter stores the conserved quantity rate generated or annihilate per volume unit. Positive values mean that the conserved quantity is created, negative ones mean that it vanishes.
|
inlineinherited |
Evaluate the source term for all phases within a given sub-control-volume.
This is the method for the case where the source term is potentially solution dependent and requires some quantities that are specific to the fully-implicit method.
element | The finite element |
fvGeometry | The finite-volume geometry |
elemVolVars | All volume variables for the element |
scv | The sub control volume |
For this method, the return parameter stores the conserved quantity rate generated or annihilate per volume unit. Positive values mean that the conserved quantity is created, negative ones mean that it vanishes. E.g. for the mass balance that would be a mass rate in \( [ kg / (m^3 \cdot s)] \).
|
inlineinherited |
Evaluate the source term for all phases within a given sub-control-volume.
globalPos | The position of the center of the finite volume for which the source term ought to be specified in global coordinates |
For this method, the values parameter stores the conserved quantity rate generated or annihilate per volume unit. Positive values mean that the conserved quantity is created, negative ones mean that it vanishes. E.g. for the mass balance that would be a mass rate in \( [ kg / (m^3 \cdot s)] \).
As a default, i.e. if the user's problem does not overload any source method return 0.0 (no source terms)
|
inlineinherited |
Return a reference to the underlying spatial parameters.
|
inlineinherited |
Return a reference to the underlying spatial parameters.