28 #ifndef EWOMS_FLASH_INTENSIVE_QUANTITIES_HH
29 #define EWOMS_FLASH_INTENSIVE_QUANTITIES_HH
37 #include <opm/material/fluidstates/CompositionalFluidState.hpp>
38 #include <opm/common/Valgrind.hpp>
40 #include <dune/common/fvector.hh>
41 #include <dune/common/fmatrix.hh>
51 template <
class TypeTag>
56 ,
public GET_PROP_TYPE(TypeTag, FluxModule)::FluxIntensiveQuantities
58 typedef typename GET_PROP_TYPE(TypeTag, DiscIntensiveQuantities) ParentType;
60 typedef typename GET_PROP_TYPE(TypeTag, ElementContext) ElementContext;
61 typedef typename GET_PROP_TYPE(TypeTag, MaterialLaw) MaterialLaw;
62 typedef typename GET_PROP_TYPE(TypeTag, MaterialLawParams) MaterialLawParams;
63 typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
64 typedef typename GET_PROP_TYPE(TypeTag, FluxModule) FluxModule;
65 typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
69 enum { cTot0Idx = Indices::cTot0Idx };
72 enum { enableDiffusion =
GET_PROP_VALUE(TypeTag, EnableDiffusion) };
74 enum { dimWorld = GridView::dimensionworld };
76 typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
77 typedef typename GET_PROP_TYPE(TypeTag, Evaluation) Evaluation;
78 typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
79 typedef typename GET_PROP_TYPE(TypeTag, FlashSolver) FlashSolver;
81 typedef Dune::FieldVector<Evaluation, numComponents> ComponentVector;
82 typedef Dune::FieldMatrix<Scalar, dimWorld, dimWorld> DimMatrix;
84 typedef typename FluxModule::FluxIntensiveQuantities FluxIntensiveQuantities;
90 typedef Opm::CompositionalFluidState<Evaluation, FluidSystem, enableEnergy>
FluidState;
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 =
EWOMS_GET_PARAM(TypeTag, Scalar, 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 Evaluation relativePermeability_[numPhases];
212 Evaluation mobility_[numPhases];
const Evaluation & mobility(unsigned phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: flashintensivequantities.hh:196
Simplifies multi-threaded capabilities.
Definition: threadmanager.hh:52
Provides the volumetric quantities required for the calculation of molecular diffusive fluxes...
Definition: diffusionmodule.hh:146
#define GET_PROP_VALUE(TypeTag, PropTagName)
Access the value attribute of a property for a type tag.
Definition: propertysystem.hh:469
Contains the classes required to consider energy as a conservation quantity in a multi-phase module...
#define GET_PROP_TYPE(TypeTag, PropTagName)
Access the type attribute of a property for a type tag.
Definition: propertysystem.hh:486
Contains the intensive quantities of the flash-based compositional multi-phase model.
Definition: flashintensivequantities.hh:52
Declares the properties required by the compositional multi-phase model based on flash calculations...
Opm::CompositionalFluidState< Evaluation, FluidSystem, enableEnergy > FluidState
The type of the object returned by the fluidState() method.
Definition: flashintensivequantities.hh:90
Provides the volumetric quantities required for the energy equation.
Definition: energymodule.hh:539
#define EWOMS_GET_PARAM(TypeTag, ParamType, ParamName)
Retrieve a runtime parameter.
Definition: parametersystem.hh:99
const Evaluation & relativePermeability(unsigned phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: flashintensivequantities.hh:190
const Evaluation & porosity() const
Returns the average porosity within the control volume.
Definition: flashintensivequantities.hh:204
void update(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Definition: flashintensivequantities.hh:102
Defines the primary variable and equation indices for the compositional multi-phase model based on fl...
Classes required for molecular diffusion.
const FluidState & fluidState() const
Returns the phase state for the control-volume.
Definition: flashintensivequantities.hh:178
const DimMatrix & intrinsicPermeability() const
Returns the intrinsic permeability tensor a degree of freedom.
Definition: flashintensivequantities.hh:184