28#ifndef EWOMS_NCP_INTENSIVE_QUANTITIES_HH
29#define EWOMS_NCP_INTENSIVE_QUANTITIES_HH
31#include <dune/common/fvector.hh>
32#include <dune/common/fmatrix.hh>
34#include <opm/material/common/Valgrind.hpp>
35#include <opm/material/constraintsolvers/CompositionFromFugacities.hpp>
36#include <opm/material/constraintsolvers/NcpFlash.hpp>
37#include <opm/material/fluidstates/CompositionalFluidState.hpp>
55template <
class TypeTag>
57 :
public GetPropType<TypeTag, Properties::DiscIntensiveQuantities>
60 ,
public GetPropType<TypeTag, Properties::FluxModule>::FluxIntensiveQuantities
75 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
76 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
77 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
78 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
79 enum { fugacity0Idx = Indices::fugacity0Idx };
80 enum { saturation0Idx = Indices::saturation0Idx };
81 enum { pressure0Idx = Indices::pressure0Idx };
82 enum { dimWorld = GridView::dimensionworld };
84 using CompositionFromFugacitiesSolver = ::Opm::CompositionFromFugacities<Scalar, FluidSystem, Evaluation>;
85 using FluidState = CompositionalFluidState<Evaluation, FluidSystem, enableEnergy>;
86 using ComponentVector = Dune::FieldVector<Evaluation, numComponents>;
87 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
90 using FluxIntensiveQuantities =
typename FluxModule::FluxIntensiveQuantities;
102 void update(
const ElementContext& elemCtx,
106 ParentType::update(elemCtx, dofIdx, timeIdx);
107 ParentType::checkDefined();
109 typename FluidSystem::template ParameterCache<Evaluation> paramCache;
110 const auto& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
113 Evaluation sumSat = 0;
114 for (
unsigned phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx) {
115 const Evaluation& val = priVars.makeEvaluation(saturation0Idx + phaseIdx, timeIdx);
116 fluidState_.setSaturation(phaseIdx, val);
119 fluidState_.setSaturation(numPhases - 1, 1.0 - sumSat);
120 Opm::Valgrind::CheckDefined(sumSat);
123 EnergyIntensiveQuantities::updateTemperatures_(fluidState_, elemCtx, dofIdx, timeIdx);
126 const auto& problem = elemCtx.problem();
127 const MaterialLawParams& materialParams =
128 problem.materialLawParams(elemCtx, dofIdx, timeIdx);
131 std::array<Evaluation, numPhases> capPress;
132 MaterialLaw::capillaryPressures(capPress, materialParams, fluidState_);
135 const Evaluation& pressure0 = priVars.makeEvaluation(pressure0Idx, timeIdx);
136 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
137 fluidState_.setPressure(phaseIdx, pressure0 + (capPress[phaseIdx] - capPress[0]));
142 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
143 fug[compIdx] = priVars.makeEvaluation(fugacity0Idx + compIdx, timeIdx);
147 const auto* hint = elemCtx.thermodynamicHint(dofIdx, timeIdx);
148 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
151 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
153 const Evaluation& moleFracIJ = hint->fluidState().moleFraction(phaseIdx, compIdx);
154 fluidState_.setMoleFraction(phaseIdx, compIdx, moleFracIJ);
158 CompositionFromFugacitiesSolver::guessInitial(fluidState_, phaseIdx, fug);
163 CompositionFromFugacitiesSolver::solve(fluidState_, paramCache, phaseIdx, fug);
167 porosity_ = problem.porosity(elemCtx, dofIdx, timeIdx);
168 Opm::Valgrind::CheckDefined(porosity_);
171 MaterialLaw::relativePermeabilities(relativePermeability_, materialParams, fluidState_);
174 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
176 const Evaluation& mu = FluidSystem::viscosity(fluidState_, paramCache, phaseIdx);
177 fluidState_.setViscosity(phaseIdx, mu);
179 mobility_[phaseIdx] = relativePermeability_[phaseIdx]/mu;
183 intrinsicPerm_ = problem.intrinsicPermeability(elemCtx, dofIdx, timeIdx);
186 FluxIntensiveQuantities::update_(elemCtx, dofIdx, timeIdx);
189 EnergyIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
192 DiffusionIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
201 {
return fluidState_; }
207 {
return intrinsicPerm_; }
213 {
return relativePermeability_[phaseIdx]; }
218 const Evaluation&
mobility(
unsigned phaseIdx)
const
219 {
return mobility_[phaseIdx]; }
225 {
return porosity_; }
232#if !defined NDEBUG && HAVE_VALGRIND
233 ParentType::checkDefined();
235 Valgrind::CheckDefined(porosity_);
236 Valgrind::CheckDefined(relativePermeability_);
238 fluidState_.checkDefined();
243 DimMatrix intrinsicPerm_;
244 FluidState fluidState_;
245 Evaluation porosity_;
246 std::array<Evaluation, numPhases> relativePermeability_{};
247 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 quantities which are are constant within a finite volume in the compositional multi-phas...
Definition: ncpintensivequantities.hh:61
void checkDefined() const
IntensiveQuantities::checkDefined.
Definition: ncpintensivequantities.hh:230
const FluidState & fluidState() const
ImmiscibleIntensiveQuantities::fluidState.
Definition: ncpintensivequantities.hh:200
const Evaluation & mobility(unsigned phaseIdx) const
ImmiscibleIntensiveQuantities::mobility.
Definition: ncpintensivequantities.hh:218
NcpIntensiveQuantities()=default
NcpIntensiveQuantities & operator=(const NcpIntensiveQuantities &other)=default
NcpIntensiveQuantities(const NcpIntensiveQuantities &other)=default
const Evaluation & porosity() const
ImmiscibleIntensiveQuantities::porosity.
Definition: ncpintensivequantities.hh:224
const Evaluation & relativePermeability(unsigned phaseIdx) const
ImmiscibleIntensiveQuantities::relativePermeability.
Definition: ncpintensivequantities.hh:212
void update(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
IntensiveQuantities::update.
Definition: ncpintensivequantities.hh:102
const DimMatrix & intrinsicPermeability() const
ImmiscibleIntensiveQuantities::intrinsicPermeability.
Definition: ncpintensivequantities.hh:206
Classes required for molecular diffusion.
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
Declare the properties used by the infrastructure code of the finite volume discretizations.
Defines the common properties required by the porous medium multi-phase models.
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