28#ifndef EWOMS_RICHARDS_INTENSIVE_QUANTITIES_HH
29#define EWOMS_RICHARDS_INTENSIVE_QUANTITIES_HH
31#include <dune/common/fmatrix.hh>
32#include <dune/common/fvector.hh>
34#include <opm/material/common/MathToolbox.hpp>
35#include <opm/material/fluidstates/ImmiscibleFluidState.hpp>
50template <
class TypeTag>
52 :
public GetPropType<TypeTag, Properties::DiscIntensiveQuantities>
53 ,
public GetPropType<TypeTag, Properties::FluxModule>::FluxIntensiveQuantities
65 enum { pressureWIdx = Indices::pressureWIdx };
66 enum { numPhases = FluidSystem::numPhases };
67 enum { liquidPhaseIdx = getPropValue<TypeTag, Properties::LiquidPhaseIndex>() };
68 enum { gasPhaseIdx = getPropValue<TypeTag, Properties::GasPhaseIndex>() };
69 enum { dimWorld = GridView::dimensionworld };
71 using FluxIntensiveQuantities =
typename FluxModule::FluxIntensiveQuantities;
72 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
73 using ScalarPhaseVector = Dune::FieldVector<Scalar, numPhases>;
74 using PhaseVector = Dune::FieldVector<Evaluation, numPhases>;
75 using Toolbox = MathToolbox<Evaluation>;
79 using FluidState = ImmiscibleFluidState<Evaluation, FluidSystem>;
90 void update(
const ElementContext& elemCtx,
unsigned dofIdx,
unsigned timeIdx)
92 ParentType::update(elemCtx, dofIdx, timeIdx);
94 const auto& T = elemCtx.problem().temperature(elemCtx, dofIdx, timeIdx);
95 fluidState_.setTemperature(T);
98 const auto& problem = elemCtx.problem();
99 const typename MaterialLaw::Params& materialParams =
100 problem.materialLawParams(elemCtx, dofIdx, timeIdx);
101 const auto& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
108 fluidState_.setSaturation(liquidPhaseIdx, 1.0);
109 fluidState_.setSaturation(gasPhaseIdx, 0.0);
110 ScalarPhaseVector pC;
111 MaterialLaw::capillaryPressures(pC, materialParams, fluidState_);
116 const Evaluation& pW = priVars.makeEvaluation(pressureWIdx, timeIdx);
117 const Evaluation pN =
118 Toolbox::max(elemCtx.problem().referencePressure(elemCtx, dofIdx, 0),
119 pW + (pC[gasPhaseIdx] - pC[liquidPhaseIdx]));
124 fluidState_.setPressure(liquidPhaseIdx, pW);
125 fluidState_.setPressure(gasPhaseIdx, pN);
128 MaterialLaw::saturations(sat, materialParams, fluidState_);
129 fluidState_.setSaturation(liquidPhaseIdx, sat[liquidPhaseIdx]);
130 fluidState_.setSaturation(gasPhaseIdx, sat[gasPhaseIdx]);
132 typename FluidSystem::template ParameterCache<Evaluation> paramCache;
133 paramCache.updateAll(fluidState_);
136 const Evaluation& mu = FluidSystem::viscosity(fluidState_, paramCache, liquidPhaseIdx);
137 fluidState_.setViscosity(liquidPhaseIdx, mu);
138 fluidState_.setViscosity(gasPhaseIdx, 1e-20);
141 const Evaluation& rho = FluidSystem::density(fluidState_, paramCache, liquidPhaseIdx);
142 fluidState_.setDensity(liquidPhaseIdx, rho);
143 fluidState_.setDensity(gasPhaseIdx, 1e-20);
146 MaterialLaw::relativePermeabilities(relativePermeability_, materialParams, fluidState_);
149 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
150 mobility_[phaseIdx] = relativePermeability_[phaseIdx] / fluidState_.viscosity(phaseIdx);
154 porosity_ = problem.porosity(elemCtx, dofIdx, timeIdx);
157 intrinsicPerm_ = problem.intrinsicPermeability(elemCtx, dofIdx, timeIdx);
160 FluxIntensiveQuantities::update_(elemCtx, dofIdx, timeIdx);
167 {
return fluidState_; }
173 {
return porosity_; }
179 {
return intrinsicPerm_; }
185 {
return relativePermeability_[phaseIdx]; }
190 const Evaluation&
mobility(
unsigned phaseIdx)
const
191 {
return mobility_[phaseIdx]; }
195 DimMatrix intrinsicPerm_;
196 std::array<Evaluation,numPhases> relativePermeability_;
197 std::array<Evaluation,numPhases> mobility_;
198 Evaluation porosity_;
Intensive quantities required by the Richards model.
Definition: richardsintensivequantities.hh:54
const Evaluation & mobility(unsigned phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: richardsintensivequantities.hh:190
const FluidState & fluidState() const
Returns the phase state for the control-volume.
Definition: richardsintensivequantities.hh:166
const Evaluation & porosity() const
Returns the average porosity within the control volume.
Definition: richardsintensivequantities.hh:172
const DimMatrix & intrinsicPermeability() const
Returns the intrinsic permeability tensor a degree of freedom.
Definition: richardsintensivequantities.hh:178
void update(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Definition: richardsintensivequantities.hh:90
RichardsIntensiveQuantities(const RichardsIntensiveQuantities &other)=default
RichardsIntensiveQuantities & operator=(const RichardsIntensiveQuantities &other)=default
ImmiscibleFluidState< Evaluation, FluidSystem > FluidState
The type returned by the fluidState() method.
Definition: richardsintensivequantities.hh:79
const Evaluation & relativePermeability(unsigned phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: richardsintensivequantities.hh:184
RichardsIntensiveQuantities()=default
Declare the properties used by the infrastructure code of the finite volume discretizations.
Definition: blackoilboundaryratevector.hh:39
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:233
Contains the property declarations for the Richards model.