Element-wise calculation of the Jacobian matrix for problems using the Richards fully implicit models. More...
#include <dumux/porousmediumflow/richards/localresidual.hh>
Element-wise calculation of the Jacobian matrix for problems using the Richards fully implicit models.
Public Member Functions | |
NumEqVector | computeStorage (const Problem &problem, const SubControlVolume &scv, const VolumeVariables &volVars) const |
Evaluates 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 |
Evaluates the mass flux over a face of a sub control volume. More... | |
template<class PartialDerivativeMatrix > | |
void | addStorageDerivatives (PartialDerivativeMatrix &partialDerivatives, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const VolumeVariables &curVolVars, const SubControlVolume &scv) const |
Adds the storage derivative. 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 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... | |
|
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.
derivativeMatrices | The matrices containing the derivatives |
problem | The problem |
element | The element |
fvGeometry | The finite volume element geometry |
curElemVolVars | The current element volume variables |
elemFluxVarsCache | The element flux variables cache |
scvf | The sub control volume face |
|
inline |
Adds flux derivatives for box method.
A | The Jacobian Matrix |
problem | The problem |
element | The element |
fvGeometry | The finite volume element geometry |
curElemVolVars | The current element volume variables |
elemFluxVarsCache | The element flux variables cache |
scvf | The sub control volume face |
|
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.
derivativeMatrices | The partial derivatives |
problem | The problem |
element | The element |
fvGeometry | The finite volume element geometry |
curElemVolVars | The current element volume variables |
elemFluxVarsCache | The element flux variables cache |
scvf | The sub control volume face |
|
inline |
Adds Robin flux derivatives for wetting and nonwetting phase.
derivativeMatrices | The matrices containing the derivatives |
problem | The problem |
element | The element |
fvGeometry | The finite volume element geometry |
curElemVolVars | The current element volume variables |
elemFluxVarsCache | The element flux variables cache |
scvf | The sub control volume face |
|
inline |
Adds source derivatives for wetting and nonwetting phase.
partialDerivatives | The partial derivatives |
problem | The problem |
element | The element |
fvGeometry | The finite volume element geometry |
curVolVars | The current volume variables |
scv | The sub control volume |
|
inline |
Adds the storage derivative.
partialDerivatives | The partial derivatives |
problem | The problem |
element | The element |
fvGeometry | The finite volume element geometry |
curVolVars | The current volume variables |
scv | The sub control volume |
|
inline |
Evaluates the mass flux over a face of a sub control volume.
problem | The problem |
element | The current element. |
fvGeometry | The finite-volume geometry |
elemVolVars | The volume variables of the current element |
scvf | The sub control volume face to compute the flux on |
elemFluxVarsCache | The cache related to flux compuation |
Add advective phase energy fluxes for the water phase only. For isothermal model the contribution is zero.
Add diffusive energy fluxes. For isothermal model the contribution is zero. The effective lambda is averaged over both fluid phases and the solid phase
|
inline |
Evaluates 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.
problem | The problem |
scv | The sub control volume |
volVars | The current or previous volVars |
The energy storage in the water, air and solid phase