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>
54template <
class TypeTag>
56 :
public GetPropType<TypeTag, Properties::DiscIntensiveQuantities>
59 ,
public GetPropType<TypeTag, Properties::FluxModule>::FluxIntensiveQuantities
72 enum { cTot0Idx = Indices::cTot0Idx };
73 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
74 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
75 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
76 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
77 enum { dimWorld = GridView::dimensionworld };
84 using ComponentVector = Dune::FieldVector<Evaluation, numComponents>;
85 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
87 using FluxIntensiveQuantities =
typename FluxModule::FluxIntensiveQuantities;
93 using FluidState = CompositionalFluidState<Evaluation, FluidSystem, enableEnergy>;
104 void update(
const ElementContext& elemCtx,
unsigned dofIdx,
unsigned timeIdx)
106 ParentType::update(elemCtx, dofIdx, timeIdx);
107 EnergyIntensiveQuantities::updateTemperatures_(fluidState_, elemCtx, dofIdx, timeIdx);
109 const auto& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
110 const auto& problem = elemCtx.problem();
111 const Scalar flashTolerance = Parameters::Get<Parameters::FlashTolerance<Scalar>>();
114 ComponentVector cTotal;
115 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
116 cTotal[compIdx] = priVars.makeEvaluation(cTot0Idx + compIdx, timeIdx);
119 const auto* hint = elemCtx.thermodynamicHint(dofIdx, timeIdx);
124 const Evaluation T = fluidState_.temperature(0);
125 fluidState_.assign(hint->fluidState());
126 fluidState_.setTemperature(T);
129 FlashSolver::guessInitial(fluidState_, cTotal);
133 typename FluidSystem::template ParameterCache<Evaluation> paramCache;
134 const MaterialLawParams& materialParams =
135 problem.materialLawParams(elemCtx, dofIdx, timeIdx);
136 FlashSolver::template solve<MaterialLaw>(fluidState_,
143 MaterialLaw::relativePermeabilities(relativePermeability_,
144 materialParams, fluidState_);
145 Valgrind::CheckDefined(relativePermeability_);
148 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
149 paramCache.updatePhase(fluidState_, phaseIdx);
151 const Evaluation& mu = FluidSystem::viscosity(fluidState_, paramCache, phaseIdx);
152 fluidState_.setViscosity(phaseIdx, mu);
154 mobility_[phaseIdx] = relativePermeability_[phaseIdx] / mu;
155 Valgrind::CheckDefined(mobility_[phaseIdx]);
163 porosity_ = problem.porosity(elemCtx, dofIdx, timeIdx);
164 Valgrind::CheckDefined(porosity_);
167 intrinsicPerm_ = problem.intrinsicPermeability(elemCtx, dofIdx, timeIdx);
170 FluxIntensiveQuantities::update_(elemCtx, dofIdx, timeIdx);
173 EnergyIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
176 DiffusionIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
183 {
return fluidState_; }
189 {
return intrinsicPerm_; }
195 {
return relativePermeability_[phaseIdx]; }
200 const Evaluation&
mobility(
unsigned phaseIdx)
const
201 {
return mobility_[phaseIdx]; }
207 {
return porosity_; }
210 DimMatrix intrinsicPerm_;
212 Evaluation porosity_;
213 std::array<Evaluation,numPhases> relativePermeability_;
214 std::array<Evaluation,numPhases> mobility_;
Provides the volumetric quantities required for the calculation of molecular diffusive fluxes.
Definition: diffusionmodule.hh:145
Provides the volumetric quantities required for the energy equation.
Definition: energymodule.hh:540
Contains the intensive quantities of the flash-based compositional multi-phase model.
Definition: flash/flashintensivequantities.hh:60
const Evaluation & mobility(unsigned phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: flash/flashintensivequantities.hh:200
const DimMatrix & intrinsicPermeability() const
Returns the intrinsic permeability tensor a degree of freedom.
Definition: flash/flashintensivequantities.hh:188
FlashIntensiveQuantities(const FlashIntensiveQuantities &other)=default
const FluidState & fluidState() const
Returns the phase state for the control-volume.
Definition: flash/flashintensivequantities.hh:182
const Evaluation & porosity() const
Returns the average porosity within the control volume.
Definition: flash/flashintensivequantities.hh:206
const Evaluation & relativePermeability(unsigned phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: flash/flashintensivequantities.hh:194
CompositionalFluidState< Evaluation, FluidSystem, enableEnergy > FluidState
The type of the object returned by the fluidState() method.
Definition: flash/flashintensivequantities.hh:93
FlashIntensiveQuantities()=default
void update(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Definition: flash/flashintensivequantities.hh:104
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: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