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

Compositional transport step in a Finite Volume discretization. More...

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

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

Description

template<class TypeTag>
class Dumux::FVTransport2P2C< 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 whole flux contribution for each cell is subdivided into a storage term, a flux term and a source term. Corresponding functions (getFlux() and getFluxOnBoundary()) are provided, internal sources are directly treated.

Template Parameters
TypeTagThe Type Tag

Classes

struct  LocalTimesteppingData
 

Public Types

using EntryType = Dune::FieldVector< Scalar, 2 >
 
using TimeStepFluxType = Dune::FieldVector< Scalar, 2 >
 

Public Member Functions

virtual void update (const Scalar t, Scalar &dt, TransportSolutionType &updateVec, bool impet=false)
 Calculate the update vector and determine timestep size This method calculates the update vector \( u \) of the discretized equation. More...
 
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 &timestepFlux, 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 &timestepFlux, 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...
 

Protected Member Functions

Problem & problem ()
 Acess function for the current problem. More...
 
void innerUpdate (TransportSolutionType &updateVec)
 

Access functions for protected variables <br>

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< LocalTimesteppingDatatimeStepData_
 Stores data for sub-time-stepping. More...
 
Scalar subCFLFactor_
 
bool localTimeStepping_
 
int verbosity_
 
static const int pressureType = getPropValue<TypeTag, Properties::PressureFormulation>()
 gives kind of pressure used ( \( 0 = p_w \), \( 1 = p_n \), \( 2 = p_{global} \)) More...
 
void updatedTargetDt_ (Scalar &dt)
 
void resetTimeStepData_ ()
 
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...
 
 FVTransport2P2C (Problem &problem)
 Constructs a FVTransport2P2C object Currently, the compositional transport scheme can not be applied with a global pressure / total velocity formulation. More...
 
virtual ~FVTransport2P2C ()
 

Member Typedef Documentation

◆ EntryType

template<class TypeTag >
using Dumux::FVTransport2P2C< TypeTag >::EntryType = Dune::FieldVector<Scalar, 2>

Type of the vector of entries

Contains the return values of the get*()-functions (matrix or right-hand side entry).

◆ TimeStepFluxType

template<class TypeTag >
using Dumux::FVTransport2P2C< TypeTag >::TimeStepFluxType = Dune::FieldVector<Scalar, 2>

Constructor & Destructor Documentation

◆ FVTransport2P2C()

template<class TypeTag >
Dumux::FVTransport2P2C< TypeTag >::FVTransport2P2C ( Problem &  problem)
inline

Constructs a FVTransport2P2C object Currently, the compositional transport scheme can not be applied with a global pressure / total velocity formulation.

Parameters
problema problem class object

◆ ~FVTransport2P2C()

template<class TypeTag >
virtual Dumux::FVTransport2P2C< TypeTag >::~FVTransport2P2C ( )
inlinevirtual

Member Function Documentation

◆ addOutputVtkFields()

template<class TypeTag >
template<class MultiWriter >
void Dumux::FVTransport2P2C< TypeTag >::addOutputVtkFields ( MultiWriter &  writer)
inline

Write transport variables into the output files.

Parameters
writerapplied VTK-writer

◆ deserializeEntity()

template<class TypeTag >
void Dumux::FVTransport2P2C< TypeTag >::deserializeEntity ( std::istream &  instream,
const Element &  element 
)
inline

Function needed for restart option of the transport model: Read in.

◆ enableLocalTimeStepping()

template<class TypeTag >
bool Dumux::FVTransport2P2C< TypeTag >::enableLocalTimeStepping ( )
inline

Function to check if local time stepping is activated.

◆ evalBoundary()

template<class TypeTag >
void Dumux::FVTransport2P2C< TypeTag >::evalBoundary ( GlobalPosition  globalPosFace,
const Intersection &  intersection,
FluidState &  BCfluidState,
PhaseVector &  pressBound 
)

Evaluate the boundary conditions.

As the transport primary variable in this formulation is the total component concentration, \( C^{\kappa} \) it seems natural that the boundary values are also total concentrations. However, as for the initial conditions, it is possible to define boundaries by means of a saturation. This choice determines which version of flash calculation is necessary to get to the composition at the boundary.

Parameters
globalPosFaceThe global position of the face of the current boundary cell
intersectionThe current intersection
BCfluidStateFluidState object that is used for the boundary
pressBoundBoundary values of phase pressures

◆ getFlux()

template<class TypeTag >
void Dumux::FVTransport2P2C< TypeTag >::getFlux ( ComponentVector &  fluxEntries,
EntryType timestepFlux,
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.

\[ - 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
fluxEntriesThe flux entries, mass influx from cell \(j\) to \(i\).
timestepFluxflow velocities for timestep estimation
intersectionThe intersection
cellDataIThe cell data for cell \(i\)

◆ getFluxOnBoundary()

template<class TypeTag >
void Dumux::FVTransport2P2C< TypeTag >::getFluxOnBoundary ( ComponentVector &  fluxEntries,
EntryType timestepFlux,
const Intersection &  intersection,
const CellData &  cellDataI 
)

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
fluxEntriesThe flux entries, mass influx from cell \(j\) to \(i\).
timestepFluxflow velocities for timestep estimation
intersectionThe intersection
cellDataIThe cell data for cell \(i\)

◆ getSource()

template<class TypeTag >
void Dumux::FVTransport2P2C< TypeTag >::getSource ( Scalar &  update,
const Element &  element,
CellData &  cellDataI 
)
inline

◆ getTransportedQuantity()

template<class TypeTag >
void Dumux::FVTransport2P2C< TypeTag >::getTransportedQuantity ( TransportSolutionType &  transportedQuantity)
inline

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.

Parameters
transportedQuantityVector of both transported components

◆ initialize()

template<class TypeTag >
void Dumux::FVTransport2P2C< TypeTag >::initialize ( )
inline

Set the initial values before the first pressure equation This method is called before first pressure equation is solved from IMPET.

◆ innerUpdate()

template<class TypeTag >
void Dumux::FVTransport2P2C< TypeTag >::innerUpdate ( TransportSolutionType &  updateVec)
protected

◆ inPhysicalRange()

template<class TypeTag >
template<class DataEntry >
bool Dumux::FVTransport2P2C< TypeTag >::inPhysicalRange ( DataEntry &  entry)
inline

Function to control the abort of the transport-sub-time-stepping depending on a physical parameter range.

Parameters
entryCell entries of the update vector

◆ problem()

template<class TypeTag >
Problem & Dumux::FVTransport2P2C< TypeTag >::problem ( )
inlineprotected

Acess function for the current problem.

◆ resetTimeStepData_()

template<class TypeTag >
void Dumux::FVTransport2P2C< TypeTag >::resetTimeStepData_ ( )
inlineprotected

◆ serializeEntity()

template<class TypeTag >
void Dumux::FVTransport2P2C< TypeTag >::serializeEntity ( std::ostream &  outstream,
const Element &  element 
)
inline

Function needed for restart option of the transport model: Write out.

◆ totalConcentration()

template<class TypeTag >
Scalar & Dumux::FVTransport2P2C< TypeTag >::totalConcentration ( int  compIdx,
int  eIdxGlobal 
)
inline

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.

Parameters
compIdxThe index of the component
eIdxGlobalThe global index of the current cell.

◆ update()

template<class TypeTag >
void Dumux::FVTransport2P2C< TypeTag >::update ( const Scalar  t,
Scalar &  dt,
TransportSolutionType &  updateVec,
bool  impet = false 
)
virtual

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
tCurrent simulation time \(\mathrm{[s]}\)
dtTime step size \(\mathrm{[s]}\)
updateVecUpdate vector, or update estimate for secants, resp. Here in \(\mathrm{[kg/m^3]}\)
impetFlag that determines if it is a real impet step or an update estimate for volume derivatives

Reimplemented in Dumux::FV2dTransport2P2CAdaptive< TypeTag >, Dumux::FV3dTransport2P2CAdaptive< TypeTag >, and Dumux::FVTransport2P2CMultiPhysics< TypeTag >.

◆ updateConcentrations()

template<class TypeTag >
void Dumux::FVTransport2P2C< TypeTag >::updateConcentrations ( TransportSolutionType &  updateVector,
Scalar  dt 
)

Updates the concentrations once an update is calculated.

This method updates both, the internal transport solution vector and the entries in the cellData.

Parameters
updateVectorUpdate vector, or update estimate for secants, resp. Here in \(\mathrm{[kg/m^3]}\)
dtTime step size \(\mathrm{[s]}\)

◆ updatedTargetDt_()

template<class TypeTag >
void Dumux::FVTransport2P2C< TypeTag >::updatedTargetDt_ ( Scalar &  dt)
protected

◆ updateTransportedQuantity() [1/2]

template<class TypeTag >
void Dumux::FVTransport2P2C< TypeTag >::updateTransportedQuantity ( TransportSolutionType &  updateVector)

Updates the transported quantity once an update is calculated.

This method updates both, the internal transport solution vector and the entries in the cellData.

Parameters
updateVectorUpdate vector, or update estimate for secants, resp. Here in \(\mathrm{[kg/m^3]}\)

◆ updateTransportedQuantity() [2/2]

template<class TypeTag >
void Dumux::FVTransport2P2C< TypeTag >::updateTransportedQuantity ( TransportSolutionType &  updateVector,
Scalar  dt 
)

Updates the transported quantity once an update is calculated.

This method updates both, the internal transport solution vector and the entries in the cellData.

Parameters
updateVectorUpdate vector, or update estimate for secants, resp. Here in \(\mathrm{[kg/m^3]}\)
dtTime step size \(\mathrm{[s]}\)

Member Data Documentation

◆ accumulatedDt_

template<class TypeTag >
Scalar Dumux::FVTransport2P2C< TypeTag >::accumulatedDt_
protected

Current time-interval in sub-time-stepping routine.

◆ averagedFaces_

template<class TypeTag >
int Dumux::FVTransport2P2C< TypeTag >::averagedFaces_
protected

number of faces were flux was restricted

◆ dtThreshold_

template<class TypeTag >
const Scalar Dumux::FVTransport2P2C< TypeTag >::dtThreshold_
protected

Threshold for sub-time-stepping routine.

◆ impet_

template<class TypeTag >
bool Dumux::FVTransport2P2C< TypeTag >::impet_
protected

indicating if we are in an estimate (false) or real impet (true) step.

◆ localTimeStepping_

template<class TypeTag >
bool Dumux::FVTransport2P2C< TypeTag >::localTimeStepping_
protected

◆ minimalBoundaryPermeability

template<class TypeTag >
Scalar Dumux::FVTransport2P2C< TypeTag >::minimalBoundaryPermeability
protected

Minimal limit for the boundary permeability.

◆ pressureType

template<class TypeTag >
const int Dumux::FVTransport2P2C< TypeTag >::pressureType = getPropValue<TypeTag, Properties::PressureFormulation>()
staticprotected

gives kind of pressure used ( \( 0 = p_w \), \( 1 = p_n \), \( 2 = p_{global} \))

◆ problem_

template<class TypeTag >
Problem& Dumux::FVTransport2P2C< TypeTag >::problem_
protected

◆ regulateBoundaryPermeability

template<class TypeTag >
bool Dumux::FVTransport2P2C< TypeTag >::regulateBoundaryPermeability
protected

Enables regulation of permeability in the direction of a Dirichlet Boundary Condition.

◆ restrictFluxInTransport_

template<class TypeTag >
int Dumux::FVTransport2P2C< TypeTag >::restrictFluxInTransport_
protected

Restriction of flux on new pressure field if direction reverses from the pressure equation.

◆ subCFLFactor_

template<class TypeTag >
Scalar Dumux::FVTransport2P2C< TypeTag >::subCFLFactor_
protected

◆ switchNormals

template<class TypeTag >
bool Dumux::FVTransport2P2C< TypeTag >::switchNormals
protected

◆ timeStepData_

template<class TypeTag >
std::vector<LocalTimesteppingData> Dumux::FVTransport2P2C< TypeTag >::timeStepData_
protected

Stores data for sub-time-stepping.

◆ totalConcentration_

template<class TypeTag >
TransportSolutionType Dumux::FVTransport2P2C< TypeTag >::totalConcentration_
protected

private vector of transported primary variables

◆ verbosity_

template<class TypeTag >
int Dumux::FVTransport2P2C< TypeTag >::verbosity_
protected

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