28#ifndef EWOMS_FLASH_INTENSIVE_QUANTITIES_HH
29#define EWOMS_FLASH_INTENSIVE_QUANTITIES_HH
37#include <opm/material/fluidstates/CompositionalFluidState.hpp>
38#include <opm/material/common/Valgrind.hpp>
40#include <dune/common/fvector.hh>
41#include <dune/common/fmatrix.hh>
51template <
class TypeTag>
53 :
public GetPropType<TypeTag, Properties::DiscIntensiveQuantities>
56 ,
public GetPropType<TypeTag, Properties::FluxModule>::FluxIntensiveQuantities
69 enum { cTot0Idx = Indices::cTot0Idx };
70 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
71 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
72 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
73 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
74 enum { dimWorld = GridView::dimensionworld };
81 using ComponentVector = Dune::FieldVector<Evaluation, numComponents>;
82 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
84 using FluxIntensiveQuantities =
typename FluxModule::FluxIntensiveQuantities;
90 using FluidState = Opm::CompositionalFluidState<Evaluation, FluidSystem, enableEnergy>;
102 void update(
const ElementContext& elemCtx,
unsigned dofIdx,
unsigned timeIdx)
104 ParentType::update(elemCtx, dofIdx, timeIdx);
105 EnergyIntensiveQuantities::updateTemperatures_(fluidState_, elemCtx, dofIdx, timeIdx);
107 const auto& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
108 const auto& problem = elemCtx.problem();
109 Scalar flashTolerance = Parameters::get<TypeTag, Properties::FlashTolerance>();
112 ComponentVector cTotal;
113 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
114 cTotal[compIdx] = priVars.makeEvaluation(cTot0Idx + compIdx, timeIdx);
116 const auto *hint = elemCtx.thermodynamicHint(dofIdx, timeIdx);
121 Evaluation T = fluidState_.temperature(0);
122 fluidState_.assign(hint->fluidState());
123 fluidState_.setTemperature(T);
126 FlashSolver::guessInitial(fluidState_, cTotal);
129 typename FluidSystem::template ParameterCache<Evaluation> paramCache;
130 const MaterialLawParams& materialParams =
131 problem.materialLawParams(elemCtx, dofIdx, timeIdx);
132 FlashSolver::template solve<MaterialLaw>(fluidState_,
139 MaterialLaw::relativePermeabilities(relativePermeability_,
140 materialParams, fluidState_);
141 Opm::Valgrind::CheckDefined(relativePermeability_);
144 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
145 paramCache.updatePhase(fluidState_, phaseIdx);
147 const Evaluation& mu = FluidSystem::viscosity(fluidState_, paramCache, phaseIdx);
148 fluidState_.setViscosity(phaseIdx, mu);
150 mobility_[phaseIdx] = relativePermeability_[phaseIdx] / mu;
151 Opm::Valgrind::CheckDefined(mobility_[phaseIdx]);
159 porosity_ = problem.porosity(elemCtx, dofIdx, timeIdx);
160 Opm::Valgrind::CheckDefined(porosity_);
163 intrinsicPerm_ = problem.intrinsicPermeability(elemCtx, dofIdx, timeIdx);
166 FluxIntensiveQuantities::update_(elemCtx, dofIdx, timeIdx);
169 EnergyIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
172 DiffusionIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
179 {
return fluidState_; }
185 {
return intrinsicPerm_; }
191 {
return relativePermeability_[phaseIdx]; }
196 const Evaluation&
mobility(
unsigned phaseIdx)
const
198 return mobility_[phaseIdx];
205 {
return porosity_; }
208 DimMatrix intrinsicPerm_;
210 Evaluation porosity_;
211 std::array<Evaluation,numPhases> relativePermeability_;
212 std::array<Evaluation,numPhases> mobility_;
Provides the volumetric quantities required for the calculation of molecular diffusive fluxes.
Definition: diffusionmodule.hh:140
Provides the volumetric quantities required for the energy equation.
Definition: energymodule.hh:530
Contains the intensive quantities of the flash-based compositional multi-phase model.
Definition: flash/flashintensivequantities.hh:57
Opm::CompositionalFluidState< Evaluation, FluidSystem, enableEnergy > FluidState
The type of the object returned by the fluidState() method.
Definition: flash/flashintensivequantities.hh:90
const Evaluation & mobility(unsigned phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: flash/flashintensivequantities.hh:196
const DimMatrix & intrinsicPermeability() const
Returns the intrinsic permeability tensor a degree of freedom.
Definition: flash/flashintensivequantities.hh:184
FlashIntensiveQuantities(const FlashIntensiveQuantities &other)=default
const FluidState & fluidState() const
Returns the phase state for the control-volume.
Definition: flash/flashintensivequantities.hh:178
const Evaluation & porosity() const
Returns the average porosity within the control volume.
Definition: flash/flashintensivequantities.hh:204
FlashIntensiveQuantities()
Definition: flash/flashintensivequantities.hh:92
const Evaluation & relativePermeability(unsigned phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: flash/flashintensivequantities.hh:190
void update(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Definition: flash/flashintensivequantities.hh:102
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.
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:242
Declares the properties required by the compositional multi-phase model based on flash calculations.