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

Element-wise calculation of the residual and its derivatives for a two-phase, incompressible test problem. More...

#include <dumux/porousmediumflow/2p/incompressiblelocalresidual.hh>

Inheritance diagram for Dumux::TwoPIncompressibleLocalResidual< TypeTag >:

Description

template<class TypeTag>
class Dumux::TwoPIncompressibleLocalResidual< TypeTag >

Element-wise calculation of the residual and its derivatives for a two-phase, incompressible test problem.

Public Member Functions

template<class PartialDerivativeMatrix >
void addStorageDerivatives (PartialDerivativeMatrix &partialDerivatives, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const VolumeVariables &curVolVars, const SubControlVolume &scv) const
 Adds storage derivatives for wetting and nonwetting phase. More...
 
template<class PartialDerivativeMatrix >
void addSourceDerivatives (PartialDerivativeMatrix &partialDerivatives, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const VolumeVariables &curVolVars, const SubControlVolume &scv) const
 Adds source derivatives for wetting and nonwetting phase. More...
 
template<class PartialDerivativeMatrices , class T = TypeTag>
std::enable_if_t< GetPropType< T, Properties::GridGeometry >::discMethod==DiscretizationMethods::cctpfa, void > addFluxDerivatives (PartialDerivativeMatrices &derivativeMatrices, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
 Adds flux derivatives for wetting and nonwetting phase for cell-centered FVM using TPFA. More...
 
template<class JacobianMatrix , class T = TypeTag>
std::enable_if_t< GetPropType< T, Properties::GridGeometry >::discMethod==DiscretizationMethods::box, void > addFluxDerivatives (JacobianMatrix &A, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
 Adds flux derivatives for wetting and nonwetting phase for box method. More...
 
template<class PartialDerivativeMatrices >
void addCCDirichletFluxDerivatives (PartialDerivativeMatrices &derivativeMatrices, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
 Adds cell-centered Dirichlet flux derivatives for wetting and nonwetting phase. More...
 
template<class PartialDerivativeMatrices >
void addRobinFluxDerivatives (PartialDerivativeMatrices &derivativeMatrices, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
 Adds Robin flux derivatives for wetting and nonwetting phase. More...
 
NumEqVector computeStorage (const Problem &problem, const SubControlVolume &scv, const VolumeVariables &volVars) const
 Evaluatex the rate of change of all conservation quantites (e.g. phase mass) within a sub-control volume of a finite volume element for the immiscible models. More...
 
NumEqVector computeFlux (const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache) const
 Evaluatex the mass flux over a face of a sub control volume. More...
 

Member Function Documentation

◆ addCCDirichletFluxDerivatives()

template<class TypeTag >
template<class PartialDerivativeMatrices >
void Dumux::TwoPIncompressibleLocalResidual< TypeTag >::addCCDirichletFluxDerivatives ( PartialDerivativeMatrices &  derivativeMatrices,
const Problem &  problem,
const Element &  element,
const FVElementGeometry &  fvGeometry,
const ElementVolumeVariables &  curElemVolVars,
const ElementFluxVariablesCache &  elemFluxVarsCache,
const SubControlVolumeFace &  scvf 
) const
inline

Adds cell-centered Dirichlet flux derivatives for wetting and nonwetting phase.

Compute derivatives for the wetting and the nonwetting phase flux with respect to \(p_w\) and \(S_n\).

Parameters
derivativeMatricesThe matrices containing the derivatives
problemThe problem
elementThe element
fvGeometryThe finite volume element geometry
curElemVolVarsThe current element volume variables
elemFluxVarsCacheThe element flux variables cache
scvfThe sub control volume face

◆ addFluxDerivatives() [1/2]

template<class TypeTag >
template<class JacobianMatrix , class T = TypeTag>
std::enable_if_t< GetPropType< T, Properties::GridGeometry >::discMethod==DiscretizationMethods::box, void > Dumux::TwoPIncompressibleLocalResidual< TypeTag >::addFluxDerivatives ( JacobianMatrix &  A,
const Problem &  problem,
const Element &  element,
const FVElementGeometry &  fvGeometry,
const ElementVolumeVariables &  curElemVolVars,
const ElementFluxVariablesCache &  elemFluxVarsCache,
const SubControlVolumeFace &  scvf 
) const
inline

Adds flux derivatives for wetting and nonwetting phase for box method.

Compute derivatives for the wetting and the nonwetting phase flux with respect to \(p_w\) and \(S_n\).

Parameters
AThe Jacobian Matrix
problemThe problem
elementThe element
fvGeometryThe finite volume element geometry
curElemVolVarsThe current element volume variables
elemFluxVarsCacheThe element flux variables cache
scvfThe sub control volume face

◆ addFluxDerivatives() [2/2]

template<class TypeTag >
template<class PartialDerivativeMatrices , class T = TypeTag>
std::enable_if_t< GetPropType< T, Properties::GridGeometry >::discMethod==DiscretizationMethods::cctpfa, void > Dumux::TwoPIncompressibleLocalResidual< TypeTag >::addFluxDerivatives ( PartialDerivativeMatrices &  derivativeMatrices,
const Problem &  problem,
const Element &  element,
const FVElementGeometry &  fvGeometry,
const ElementVolumeVariables &  curElemVolVars,
const ElementFluxVariablesCache &  elemFluxVarsCache,
const SubControlVolumeFace &  scvf 
) const
inline

Adds flux derivatives for wetting and nonwetting phase for cell-centered FVM using TPFA.

Compute derivatives for the wetting and the nonwetting phase flux with respect to \(p_w\) and \(S_n\).

Parameters
derivativeMatricesThe partial derivatives
problemThe problem
elementThe element
fvGeometryThe finite volume element geometry
curElemVolVarsThe current element volume variables
elemFluxVarsCacheThe element flux variables cache
scvfThe sub control volume face

◆ addRobinFluxDerivatives()

template<class TypeTag >
template<class PartialDerivativeMatrices >
void Dumux::TwoPIncompressibleLocalResidual< TypeTag >::addRobinFluxDerivatives ( PartialDerivativeMatrices &  derivativeMatrices,
const Problem &  problem,
const Element &  element,
const FVElementGeometry &  fvGeometry,
const ElementVolumeVariables &  curElemVolVars,
const ElementFluxVariablesCache &  elemFluxVarsCache,
const SubControlVolumeFace &  scvf 
) const
inline

Adds Robin flux derivatives for wetting and nonwetting phase.

Parameters
derivativeMatricesThe matrices containing the derivatives
problemThe problem
elementThe element
fvGeometryThe finite volume element geometry
curElemVolVarsThe current element volume variables
elemFluxVarsCacheThe element flux variables cache
scvfThe sub control volume face

◆ addSourceDerivatives()

template<class TypeTag >
template<class PartialDerivativeMatrix >
void Dumux::TwoPIncompressibleLocalResidual< TypeTag >::addSourceDerivatives ( PartialDerivativeMatrix &  partialDerivatives,
const Problem &  problem,
const Element &  element,
const FVElementGeometry &  fvGeometry,
const VolumeVariables &  curVolVars,
const SubControlVolume &  scv 
) const
inline

Adds source derivatives for wetting and nonwetting phase.

Parameters
partialDerivativesThe partial derivatives
problemThe problem
elementThe element
fvGeometryThe finite volume element geometry
curVolVarsThe current volume variables
scvThe sub control volume

◆ addStorageDerivatives()

template<class TypeTag >
template<class PartialDerivativeMatrix >
void Dumux::TwoPIncompressibleLocalResidual< TypeTag >::addStorageDerivatives ( PartialDerivativeMatrix &  partialDerivatives,
const Problem &  problem,
const Element &  element,
const FVElementGeometry &  fvGeometry,
const VolumeVariables &  curVolVars,
const SubControlVolume &  scv 
) const
inline

Adds storage derivatives for wetting and nonwetting phase.

Compute storage derivatives for the wetting and the nonwetting phase with respect to \(p_w\) and \(S_n\).

Parameters
partialDerivativesThe partial derivatives
problemThe problem
elementThe element
fvGeometryThe finite volume element geometry
curVolVarsThe current volume variables
scvThe sub control volume

◆ computeFlux()

template<class TypeTag >
NumEqVector Dumux::ImmiscibleLocalResidual< TypeTag >::computeFlux ( const Problem &  problem,
const Element &  element,
const FVElementGeometry &  fvGeometry,
const ElementVolumeVariables &  elemVolVars,
const SubControlVolumeFace &  scvf,
const ElementFluxVariablesCache &  elemFluxVarsCache 
) const
inlineinherited

Evaluatex the mass flux over a face of a sub control volume.

Parameters
problemThe problem
elementThe element
fvGeometryThe finite volume geometry context
elemVolVarsThe volume variables for all flux stencil elements
scvfThe sub control volume face to compute the flux on
elemFluxVarsCacheThe cache related to flux computation

Add advective phase energy fluxes. For isothermal model the contribution is zero.

Add diffusive energy fluxes. For isothermal model the contribution is zero.

◆ computeStorage()

template<class TypeTag >
NumEqVector Dumux::ImmiscibleLocalResidual< TypeTag >::computeStorage ( const Problem &  problem,
const SubControlVolume &  scv,
const VolumeVariables &  volVars 
) const
inlineinherited

Evaluatex the rate of change of all conservation quantites (e.g. phase mass) within a sub-control volume of a finite volume element for the immiscible models.

Parameters
problemTODO docme!
scvThe sub control volume
volVarsThe current or previous volVars
Note
This function should not include the source and sink terms.
The volVars can be different to allow computing the implicit euler time derivative here

The energy storage in the fluid phase with index phaseIdx

The energy storage in the solid matrix


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