template<class TypeTag>
class Dumux::FVTransport2P2CMultiPhysics< TypeTag >
Compositional transport step in a finite volume discretization.
The finite volume model for the solution of the transport equation for compositional two-phase flow.
\[ \frac{\partial C^\kappa}{\partial t} = - \nabla \cdot \left( \sum_{\alpha} X^{\kappa}_{\alpha} \varrho_{alpha} \bf{v}_{\alpha}\right) + q^{\kappa}, \]
where \( \bf{v}_{\alpha} = - \lambda_{\alpha} \bf{K} \left(\nabla p_{\alpha} + \rho_{\alpha} \bf{g} \right) \). \( p_{\alpha} \) denotes the phase pressure, \( \bf{K} \) the absolute permeability, \( \lambda_{\alpha} \) the phase mobility, \( \rho_{\alpha} \) the phase density and \( \bf{g} \) the gravity constant and \( C^{\kappa} \) the total Component concentration.
The model domain is automatically divided into a single-phase and a two-phase domain. As the flux computation is relatively cheap, the same method is used for the real transport step independently of the subdomain. The pressure equation does not need any derivatives in simple subdomains, therefore in the transport estimate step inter-cell fluxes in the simple subdomain are omitted.
- Template Parameters
-
|
virtual void | update (const Scalar t, Scalar &dt, TransportSolutionType &updateVec, bool impet=false) |
| Calculate the update vector and determine timestep size. More...
|
|
| FVTransport2P2CMultiPhysics (Problem &problem) |
| Constructs a FVTransport2P2CMultiPhysics object. More...
|
|
virtual | ~FVTransport2P2CMultiPhysics () |
|
void | updateTransportedQuantity (TransportSolutionType &updateVector) |
| Updates the transported quantity once an update is calculated. More...
|
|
void | updateTransportedQuantity (TransportSolutionType &updateVector, Scalar dt) |
| Updates the transported quantity once an update is calculated. More...
|
|
void | updateConcentrations (TransportSolutionType &updateVector, Scalar dt) |
| Updates the concentrations once an update is calculated. More...
|
|
void | getFlux (ComponentVector &fluxEntries, EntryType ×tepFlux, const Intersection &intersection, CellData &cellDataI) |
| Get flux at an interface between two cells The flux through \( \gamma \) is calculated according to the underlying pressure field, calculated by the pressure model. More...
|
|
void | getFluxOnBoundary (ComponentVector &fluxEntries, EntryType ×tepFlux, const Intersection &intersection, const CellData &cellDataI) |
| Get flux on Boundary. More...
|
|
void | evalBoundary (GlobalPosition globalPosFace, const Intersection &intersection, FluidState &BCfluidState, PhaseVector &pressBound) |
| Evaluate the boundary conditions. More...
|
|
void | initialize () |
| Set the initial values before the first pressure equation This method is called before first pressure equation is solved from IMPET. More...
|
|
template<class MultiWriter > |
void | addOutputVtkFields (MultiWriter &writer) |
| Write transport variables into the output files. More...
|
|
void | serializeEntity (std::ostream &outstream, const Element &element) |
| Function needed for restart option of the transport model: Write out. More...
|
|
void | deserializeEntity (std::istream &instream, const Element &element) |
| Function needed for restart option of the transport model: Read in. More...
|
|
|
void | getTransportedQuantity (TransportSolutionType &transportedQuantity) |
| Return the vector of the transported quantity For an immiscible IMPES scheme, this is the saturation. For compositional simulations, however, the total concentration of all components is transported. More...
|
|
Scalar & | totalConcentration (int compIdx, int eIdxGlobal) |
| Return the the total concentration stored in the transport vector To get real cell values, do not acess this method, but rather call the respective function in the cell data object. More...
|
|
void | getSource (Scalar &update, const Element &element, CellData &cellDataI) |
|
template<class DataEntry > |
bool | inPhysicalRange (DataEntry &entry) |
| Function to control the abort of the transport-sub-time-stepping depending on a physical parameter range. More...
|
|
bool | enableLocalTimeStepping () |
| Function to check if local time stepping is activated. More...
|
|
TransportSolutionType | totalConcentration_ |
| private vector of transported primary variables More...
|
|
Problem & | problem_ |
|
bool | impet_ |
| indicating if we are in an estimate (false) or real impet (true) step. More...
|
|
int | averagedFaces_ |
| number of faces were flux was restricted More...
|
|
int | restrictFluxInTransport_ |
| Restriction of flux on new pressure field if direction reverses from the pressure equation. More...
|
|
bool | regulateBoundaryPermeability |
| Enables regulation of permeability in the direction of a Dirichlet Boundary Condition. More...
|
|
Scalar | minimalBoundaryPermeability |
| Minimal limit for the boundary permeability. More...
|
|
bool | switchNormals |
|
Scalar | accumulatedDt_ |
| Current time-interval in sub-time-stepping routine. More...
|
|
const Scalar | dtThreshold_ |
| Threshold for sub-time-stepping routine. More...
|
|
std::vector< LocalTimesteppingData > | timeStepData_ |
| Stores data for sub-time-stepping. More...
|
|
Scalar | subCFLFactor_ |
|
bool | localTimeStepping_ |
|
int | verbosity_ |
|
static const int | pressureType = GET_PROP_VALUE(TypeTag, PressureFormulation) |
| gives kind of pressure used ( \( 0 = p_w \), \( 1 = p_n \), \( 2 = p_{global} \)) More...
|
|
void | updatedTargetDt_ (Scalar &dt) |
|
void | resetTimeStepData_ () |
|
template<class TypeTag >
void Dumux::FVTransport2P2C< TypeTag >::getFlux |
( |
ComponentVector & |
fluxEntries, |
|
|
EntryType & |
timestepFlux, |
|
|
const Intersection & |
intersection, |
|
|
CellData & |
cellDataI |
|
) |
| |
|
inherited |
Get flux at an interface between two cells The flux through \( \gamma \) is calculated according to the underlying pressure field, calculated by the pressure model.
\[ - A_{\gamma} \mathbf{n}^T_{\gamma} \mathbf{K} \sum_{\alpha} \varrho_{\alpha} \lambda_{\alpha} \mathbf{d}_{ij} \left( \frac{p_{\alpha,j}^t - p^{t}_{\alpha,i}}{\Delta x} + \varrho_{\alpha} \mathbf{g}^T \mathbf{d}_{ij} \right) \sum_{\kappa} X^{\kappa}_{\alpha} \]
Here, \( \mathbf{d}_{ij} \) is the normalized vector connecting the cell centers, and \( \mathbf{n}_{\gamma} \) represents the normal of the face \( \gamma \). Due to the nature of the primay Variable, the (volume-)specific total mass concentration, this represents a mass flux per cell volume.
- Parameters
-
fluxEntries | The flux entries, mass influx from cell \(j\) to \(i\). |
timestepFlux | flow velocities for timestep estimation |
intersection | The intersection |
cellDataI | The cell data for cell \(i\) |
template<class TypeTag >
void Dumux::FVTransport2P2C< TypeTag >::getFluxOnBoundary |
( |
ComponentVector & |
fluxEntries, |
|
|
EntryType & |
timestepFlux, |
|
|
const Intersection & |
intersection, |
|
|
const CellData & |
cellDataI |
|
) |
| |
|
inherited |
Get flux on Boundary.
The flux through \( \gamma \) is calculated according to the underlying pressure field, calculated by the pressure model.
\[ - A_{\gamma} \mathbf{n}^T_{\gamma} \mathbf{K} \mathbf{d}_{i-Boundary} \sum_{\alpha} \varrho_{\alpha} \lambda_{\alpha} \sum_{\kappa} X^{\kappa}_{\alpha} \left( \frac{p_{\alpha,Boundary}^t - p^{t}_{\alpha,i}}{\Delta x} + \varrho_{\alpha}\mathbf{g}^T \mathbf{d}_{i-Boundary}\right) \]
Here, \( \mathbf{u} \) is the normalized vector connecting the cell centers, and \( \mathbf{n}_{\gamma_{ij}} \) represents the normal of the face \( \gamma{ij} \).
- Parameters
-
fluxEntries | The flux entries, mass influx from cell \(j\) to \(i\). |
timestepFlux | flow velocities for timestep estimation |
intersection | The intersection |
cellDataI | The cell data for cell \(i\) |
Calculate the update vector and determine timestep size.
This method calculates the update vector \( u \) of the discretized equation
\[ C^{\kappa , new} = C^{\kappa , old} + u, \]
where \( u = \sum_{\gamma} \boldsymbol{v}_{\alpha} * \varrho_{\alpha} * X^{\kappa}_{\alpha} * \boldsymbol{n} * A_{\gamma} \), \( \boldsymbol{n} \) is the face normal and \( A_{\gamma} \) is the face area of face \( \gamma \).
In addition to the update vector, the recommended time step size dt is calculated employing a CFL condition.
- Parameters
-
| t | Current simulation time \(\mathrm{[s]}\) |
[out] | dt | Time step size \(\mathrm{[s]}\) |
[out] | updateVec | Update vector, or update estimate for secants, resp. Here in \(\mathrm{[kg/m^3]}\) |
| impet | Flag that determines if it is a real impet step or an update estimate for volume derivatives |
Reimplemented from Dumux::FVTransport2P2C< TypeTag >.