28#ifndef EWOMS_NCP_INTENSIVE_QUANTITIES_HH
29#define EWOMS_NCP_INTENSIVE_QUANTITIES_HH
36#include <opm/material/constraintsolvers/NcpFlash.hpp>
37#include <opm/material/fluidstates/CompositionalFluidState.hpp>
38#include <opm/material/constraintsolvers/CompositionFromFugacities.hpp>
39#include <opm/material/common/Valgrind.hpp>
41#include <dune/common/fvector.hh>
42#include <dune/common/fmatrix.hh>
52template <
class TypeTag>
54 :
public GetPropType<TypeTag, Properties::DiscIntensiveQuantities>
57 ,
public GetPropType<TypeTag, Properties::FluxModule>::FluxIntensiveQuantities
72 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
73 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
74 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
75 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
76 enum { fugacity0Idx = Indices::fugacity0Idx };
77 enum { saturation0Idx = Indices::saturation0Idx };
78 enum { pressure0Idx = Indices::pressure0Idx };
79 enum { dimWorld = GridView::dimensionworld };
81 using CompositionFromFugacitiesSolver = Opm::CompositionFromFugacities<Scalar, FluidSystem, Evaluation>;
82 using FluidState = Opm::CompositionalFluidState<Evaluation, FluidSystem, enableEnergy>;
83 using ComponentVector = Dune::FieldVector<Evaluation, numComponents>;
84 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
87 using FluxIntensiveQuantities =
typename FluxModule::FluxIntensiveQuantities;
100 void update(
const ElementContext& elemCtx,
104 ParentType::update(elemCtx, dofIdx, timeIdx);
105 ParentType::checkDefined();
107 typename FluidSystem::template ParameterCache<Evaluation> paramCache;
108 const auto& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
111 Evaluation sumSat = 0;
112 for (
unsigned phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx) {
113 const Evaluation& val = priVars.makeEvaluation(saturation0Idx + phaseIdx, timeIdx);
114 fluidState_.setSaturation(phaseIdx, val);
117 fluidState_.setSaturation(numPhases - 1, 1.0 - sumSat);
118 Opm::Valgrind::CheckDefined(sumSat);
121 EnergyIntensiveQuantities::updateTemperatures_(fluidState_, elemCtx, dofIdx, timeIdx);
124 const auto& problem = elemCtx.problem();
125 const MaterialLawParams& materialParams =
126 problem.materialLawParams(elemCtx, dofIdx, timeIdx);
128 Evaluation capPress[numPhases];
129 MaterialLaw::capillaryPressures(capPress, materialParams, fluidState_);
131 const Evaluation& pressure0 = priVars.makeEvaluation(pressure0Idx, timeIdx);
132 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
133 fluidState_.setPressure(phaseIdx, pressure0 + (capPress[phaseIdx] - capPress[0]));
137 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
138 fug[compIdx] = priVars.makeEvaluation(fugacity0Idx + compIdx, timeIdx);
141 const auto *hint = elemCtx.thermodynamicHint(dofIdx, timeIdx);
142 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
145 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
147 const Evaluation& moleFracIJ = hint->fluidState().moleFraction(phaseIdx, compIdx);
148 fluidState_.setMoleFraction(phaseIdx, compIdx, moleFracIJ);
152 CompositionFromFugacitiesSolver::guessInitial(fluidState_, phaseIdx, fug);
156 CompositionFromFugacitiesSolver::solve(fluidState_, paramCache, phaseIdx, fug);
160 porosity_ = problem.porosity(elemCtx, dofIdx, timeIdx);
161 Opm::Valgrind::CheckDefined(porosity_);
164 MaterialLaw::relativePermeabilities(relativePermeability_, materialParams, fluidState_);
167 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
169 const Evaluation& mu = FluidSystem::viscosity(fluidState_, paramCache, phaseIdx);
170 fluidState_.setViscosity(phaseIdx, mu);
172 mobility_[phaseIdx] = relativePermeability_[phaseIdx]/mu;
176 intrinsicPerm_ = problem.intrinsicPermeability(elemCtx, dofIdx, timeIdx);
179 FluxIntensiveQuantities::update_(elemCtx, dofIdx, timeIdx);
182 EnergyIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
185 DiffusionIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
194 {
return fluidState_; }
200 {
return intrinsicPerm_; }
206 {
return relativePermeability_[phaseIdx]; }
211 const Evaluation&
mobility(
unsigned phaseIdx)
const
212 {
return mobility_[phaseIdx]; }
218 {
return porosity_; }
225#if !defined NDEBUG && HAVE_VALGRIND
226 ParentType::checkDefined();
228 Opm::Valgrind::CheckDefined(porosity_);
229 Opm::Valgrind::CheckDefined(relativePermeability_);
231 fluidState_.checkDefined();
236 DimMatrix intrinsicPerm_;
237 FluidState fluidState_;
238 Evaluation porosity_;
239 Evaluation relativePermeability_[numPhases];
240 Evaluation mobility_[numPhases];
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 quantities which are are constant within a finite volume in the compositional multi-phas...
Definition: ncpintensivequantities.hh:58
void checkDefined() const
IntensiveQuantities::checkDefined.
Definition: ncpintensivequantities.hh:223
const FluidState & fluidState() const
ImmiscibleIntensiveQuantities::fluidState.
Definition: ncpintensivequantities.hh:193
const Evaluation & mobility(unsigned phaseIdx) const
ImmiscibleIntensiveQuantities::mobility.
Definition: ncpintensivequantities.hh:211
NcpIntensiveQuantities & operator=(const NcpIntensiveQuantities &other)=default
NcpIntensiveQuantities(const NcpIntensiveQuantities &other)=default
NcpIntensiveQuantities()
Definition: ncpintensivequantities.hh:90
const Evaluation & porosity() const
ImmiscibleIntensiveQuantities::porosity.
Definition: ncpintensivequantities.hh:217
const Evaluation & relativePermeability(unsigned phaseIdx) const
ImmiscibleIntensiveQuantities::relativePermeability.
Definition: ncpintensivequantities.hh:205
void update(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
IntensiveQuantities::update.
Definition: ncpintensivequantities.hh:100
const DimMatrix & intrinsicPermeability() const
ImmiscibleIntensiveQuantities::intrinsicPermeability.
Definition: ncpintensivequantities.hh:199
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:235
Declares the properties required for the NCP compositional multi-phase model.