28#ifndef EWOMS_FLASH_INTENSIVE_QUANTITIES_HH
29#define EWOMS_FLASH_INTENSIVE_QUANTITIES_HH
31#include <dune/common/fmatrix.hh>
32#include <dune/common/fvector.hh>
34#include <opm/material/common/Valgrind.hpp>
35#include <opm/material/fluidstates/CompositionalFluidState.hpp>
52template <
class TypeTag>
54 :
public GetPropType<TypeTag, Properties::DiscIntensiveQuantities>
57 ,
public GetPropType<TypeTag, Properties::FluxModule>::FluxIntensiveQuantities
70 enum { cTot0Idx = Indices::cTot0Idx };
71 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
72 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
73 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
74 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
75 enum { dimWorld = GridView::dimensionworld };
82 using ComponentVector = Dune::FieldVector<Evaluation, numComponents>;
83 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
85 using FluxIntensiveQuantities =
typename FluxModule::FluxIntensiveQuantities;
91 using FluidState = Opm::CompositionalFluidState<Evaluation, FluidSystem, enableEnergy>;
103 void update(
const ElementContext& elemCtx,
unsigned dofIdx,
unsigned timeIdx)
105 ParentType::update(elemCtx, dofIdx, timeIdx);
106 EnergyIntensiveQuantities::updateTemperatures_(fluidState_, elemCtx, dofIdx, timeIdx);
108 const auto& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
109 const auto& problem = elemCtx.problem();
110 Scalar flashTolerance = Parameters::Get<Parameters::FlashTolerance<Scalar>>();
113 ComponentVector cTotal;
114 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
115 cTotal[compIdx] = priVars.makeEvaluation(cTot0Idx + compIdx, timeIdx);
117 const auto *hint = elemCtx.thermodynamicHint(dofIdx, timeIdx);
122 Evaluation T = fluidState_.temperature(0);
123 fluidState_.assign(hint->fluidState());
124 fluidState_.setTemperature(T);
127 FlashSolver::guessInitial(fluidState_, cTotal);
130 typename FluidSystem::template ParameterCache<Evaluation> paramCache;
131 const MaterialLawParams& materialParams =
132 problem.materialLawParams(elemCtx, dofIdx, timeIdx);
133 FlashSolver::template solve<MaterialLaw>(fluidState_,
140 MaterialLaw::relativePermeabilities(relativePermeability_,
141 materialParams, fluidState_);
142 Opm::Valgrind::CheckDefined(relativePermeability_);
145 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
146 paramCache.updatePhase(fluidState_, phaseIdx);
148 const Evaluation& mu = FluidSystem::viscosity(fluidState_, paramCache, phaseIdx);
149 fluidState_.setViscosity(phaseIdx, mu);
151 mobility_[phaseIdx] = relativePermeability_[phaseIdx] / mu;
152 Opm::Valgrind::CheckDefined(mobility_[phaseIdx]);
160 porosity_ = problem.porosity(elemCtx, dofIdx, timeIdx);
161 Opm::Valgrind::CheckDefined(porosity_);
164 intrinsicPerm_ = problem.intrinsicPermeability(elemCtx, dofIdx, timeIdx);
167 FluxIntensiveQuantities::update_(elemCtx, dofIdx, timeIdx);
170 EnergyIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
173 DiffusionIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
180 {
return fluidState_; }
186 {
return intrinsicPerm_; }
192 {
return relativePermeability_[phaseIdx]; }
197 const Evaluation&
mobility(
unsigned phaseIdx)
const
199 return mobility_[phaseIdx];
206 {
return porosity_; }
209 DimMatrix intrinsicPerm_;
211 Evaluation porosity_;
212 std::array<Evaluation,numPhases> relativePermeability_;
213 std::array<Evaluation,numPhases> mobility_;
Provides the volumetric quantities required for the calculation of molecular diffusive fluxes.
Definition: diffusionmodule.hh:141
Provides the volumetric quantities required for the energy equation.
Definition: energymodule.hh:532
Contains the intensive quantities of the flash-based compositional multi-phase model.
Definition: flash/flashintensivequantities.hh:58
Opm::CompositionalFluidState< Evaluation, FluidSystem, enableEnergy > FluidState
The type of the object returned by the fluidState() method.
Definition: flash/flashintensivequantities.hh:91
const Evaluation & mobility(unsigned phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: flash/flashintensivequantities.hh:197
const DimMatrix & intrinsicPermeability() const
Returns the intrinsic permeability tensor a degree of freedom.
Definition: flash/flashintensivequantities.hh:185
FlashIntensiveQuantities(const FlashIntensiveQuantities &other)=default
const FluidState & fluidState() const
Returns the phase state for the control-volume.
Definition: flash/flashintensivequantities.hh:179
const Evaluation & porosity() const
Returns the average porosity within the control volume.
Definition: flash/flashintensivequantities.hh:205
FlashIntensiveQuantities()
Definition: flash/flashintensivequantities.hh:93
const Evaluation & relativePermeability(unsigned phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: flash/flashintensivequantities.hh:191
void update(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Definition: flash/flashintensivequantities.hh:103
FlashIntensiveQuantities & operator=(const FlashIntensiveQuantities &other)=default
Classes required for molecular diffusion.
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
Declares the parameters for the compositional multi-phase model based on flash calculations.
Declares the properties required by the compositional multi-phase model based on flash calculations.
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