28#ifndef EWOMS_RICHARDS_INTENSIVE_QUANTITIES_HH
29#define EWOMS_RICHARDS_INTENSIVE_QUANTITIES_HH
33#include <opm/material/fluidstates/ImmiscibleFluidState.hpp>
35#include <dune/common/fvector.hh>
36#include <dune/common/fmatrix.hh>
46template <
class TypeTag>
48 :
public GetPropType<TypeTag, Properties::DiscIntensiveQuantities>
49 ,
public GetPropType<TypeTag, Properties::FluxModule>::FluxIntensiveQuantities
61 enum { pressureWIdx = Indices::pressureWIdx };
62 enum { numPhases = FluidSystem::numPhases };
63 enum { liquidPhaseIdx = getPropValue<TypeTag, Properties::LiquidPhaseIndex>() };
64 enum { gasPhaseIdx = getPropValue<TypeTag, Properties::GasPhaseIndex>() };
65 enum { dimWorld = GridView::dimensionworld };
67 using FluxIntensiveQuantities =
typename FluxModule::FluxIntensiveQuantities;
68 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
69 using ScalarPhaseVector = Dune::FieldVector<Scalar, numPhases>;
70 using PhaseVector = Dune::FieldVector<Evaluation, numPhases>;
71 using Toolbox = Opm::MathToolbox<Evaluation>;
75 using FluidState = Opm::ImmiscibleFluidState<Evaluation, FluidSystem>;
87 void update(
const ElementContext& elemCtx,
unsigned dofIdx,
unsigned timeIdx)
89 ParentType::update(elemCtx, dofIdx, timeIdx);
91 const auto& T = elemCtx.problem().temperature(elemCtx, dofIdx, timeIdx);
92 fluidState_.setTemperature(T);
95 const auto& problem = elemCtx.problem();
96 const typename MaterialLaw::Params& materialParams =
97 problem.materialLawParams(elemCtx, dofIdx, timeIdx);
98 const auto& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
105 fluidState_.setSaturation(liquidPhaseIdx, 1.0);
106 fluidState_.setSaturation(gasPhaseIdx, 0.0);
107 ScalarPhaseVector pC;
108 MaterialLaw::capillaryPressures(pC, materialParams, fluidState_);
113 const Evaluation& pW = priVars.makeEvaluation(pressureWIdx, timeIdx);
115 Toolbox::max(elemCtx.problem().referencePressure(elemCtx, dofIdx, 0),
116 pW + (pC[gasPhaseIdx] - pC[liquidPhaseIdx]));
121 fluidState_.setPressure(liquidPhaseIdx, pW);
122 fluidState_.setPressure(gasPhaseIdx, pN);
125 MaterialLaw::saturations(sat, materialParams, fluidState_);
126 fluidState_.setSaturation(liquidPhaseIdx, sat[liquidPhaseIdx]);
127 fluidState_.setSaturation(gasPhaseIdx, sat[gasPhaseIdx]);
129 typename FluidSystem::template ParameterCache<Evaluation> paramCache;
130 paramCache.updateAll(fluidState_);
133 const Evaluation& mu = FluidSystem::viscosity(fluidState_, paramCache, liquidPhaseIdx);
134 fluidState_.setViscosity(liquidPhaseIdx, mu);
135 fluidState_.setViscosity(gasPhaseIdx, 1e-20);
138 const Evaluation& rho = FluidSystem::density(fluidState_, paramCache, liquidPhaseIdx);
139 fluidState_.setDensity(liquidPhaseIdx, rho);
140 fluidState_.setDensity(gasPhaseIdx, 1e-20);
143 MaterialLaw::relativePermeabilities(relativePermeability_, materialParams, fluidState_);
146 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
147 mobility_[phaseIdx] = relativePermeability_[phaseIdx]/fluidState_.viscosity(phaseIdx);
150 porosity_ = problem.porosity(elemCtx, dofIdx, timeIdx);
153 intrinsicPerm_ = problem.intrinsicPermeability(elemCtx, dofIdx, timeIdx);
156 FluxIntensiveQuantities::update_(elemCtx, dofIdx, timeIdx);
163 {
return fluidState_; }
169 {
return porosity_; }
175 {
return intrinsicPerm_; }
181 {
return relativePermeability_[phaseIdx]; }
186 const Evaluation&
mobility(
unsigned phaseIdx)
const
187 {
return mobility_[phaseIdx]; }
191 DimMatrix intrinsicPerm_;
192 std::array<Evaluation,numPhases> relativePermeability_;
193 std::array<Evaluation,numPhases> mobility_;
194 Evaluation porosity_;
Intensive quantities required by the Richards model.
Definition: richardsintensivequantities.hh:50
const Evaluation & mobility(unsigned phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: richardsintensivequantities.hh:186
const FluidState & fluidState() const
Returns the phase state for the control-volume.
Definition: richardsintensivequantities.hh:162
const Evaluation & porosity() const
Returns the average porosity within the control volume.
Definition: richardsintensivequantities.hh:168
RichardsIntensiveQuantities()
Definition: richardsintensivequantities.hh:77
const DimMatrix & intrinsicPermeability() const
Returns the intrinsic permeability tensor a degree of freedom.
Definition: richardsintensivequantities.hh:174
void update(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Definition: richardsintensivequantities.hh:87
Opm::ImmiscibleFluidState< Evaluation, FluidSystem > FluidState
The type returned by the fluidState() method.
Definition: richardsintensivequantities.hh:75
RichardsIntensiveQuantities(const RichardsIntensiveQuantities &other)=default
RichardsIntensiveQuantities & operator=(const RichardsIntensiveQuantities &other)=default
const Evaluation & relativePermeability(unsigned phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: richardsintensivequantities.hh:180
Definition: blackoilboundaryratevector.hh:37
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:235
Contains the property declarations for the Richards model.