28#ifndef EWOMS_BLACK_OIL_RATE_VECTOR_HH
29#define EWOMS_BLACK_OIL_RATE_VECTOR_HH
31#include <dune/common/fvector.hh>
33#include <opm/material/common/MathToolbox.hpp>
34#include <opm/material/common/Valgrind.hpp>
35#include <opm/material/constraintsolvers/NcpFlash.hpp>
58template <
class TypeTag>
60 :
public Dune::FieldVector<GetPropType<TypeTag, Properties::Evaluation>,
61 getPropValue<TypeTag, Properties::NumEq>()>
74 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
75 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
76 enum { conti0EqIdx = Indices::conti0EqIdx };
77 enum { contiEnergyEqIdx = Indices::contiEnergyEqIdx };
78 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
79 enum { enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>() };
80 enum { enablePolymer = getPropValue<TypeTag, Properties::EnablePolymer>() };
81 enum { enablePolymerMolarWeight = getPropValue<TypeTag, Properties::EnablePolymerMW>() };
82 enum { enableFoam = getPropValue<TypeTag, Properties::EnableFoam>() };
83 enum { enableBrine = getPropValue<TypeTag, Properties::EnableBrine>() };
84 enum { enableMICP = getPropValue<TypeTag, Properties::EnableMICP>() };
85 using Toolbox = MathToolbox<Evaluation>;
86 using ParentType = Dune::FieldVector<Evaluation, numEq>;
90 { Valgrind::SetUndefined(*
this); }
101 void setMassRate(
const ParentType& value,
unsigned pvtRegionIdx = 0)
103 ParentType::operator=(value);
106 if constexpr (getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>()) {
107 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
108 (*this)[Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx)] /=
109 FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, pvtRegionIdx);
111 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
112 (*this)[Indices::canonicalToActiveComponentIndex(FluidSystem::oilCompIdx)] /=
113 FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, pvtRegionIdx);
115 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
116 (*this)[Indices::canonicalToActiveComponentIndex(FluidSystem::waterCompIdx)] /=
117 FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, pvtRegionIdx);
119 if constexpr (enableSolvent) {
121 (*this)[Indices::contiSolventEqIdx] /=
122 solventPvt.referenceDensity(pvtRegionIdx);
133 ParentType::operator=(value);
136 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
137 (*this)[conti0EqIdx + compIdx] *= FluidSystem::molarMass(compIdx, pvtRegionIdx);
141 (*this)[Indices::contiSolventEqIdx] *= solventPvt.molarMass(pvtRegionIdx);
143 if constexpr (enablePolymer) {
144 if constexpr (enablePolymerMolarWeight) {
145 throw std::logic_error(
"Set molar rate with polymer weight tracking not implemented");
151 if constexpr (enableFoam) {
152 throw std::logic_error(
"setMolarRate() not implemented for foam");
155 if constexpr (enableBrine) {
156 throw std::logic_error(
"setMolarRate() not implemented for salt water");
159 if constexpr (enableMICP) {
160 throw std::logic_error(
"setMolarRate() not implemented for MICP");
164 if constexpr (getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>()) {
165 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
166 (*this)[Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx)] /=
167 FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, pvtRegionIdx);
169 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
170 (*this)[Indices::canonicalToActiveComponentIndex(FluidSystem::oilCompIdx)] /=
171 FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, pvtRegionIdx);
173 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
174 (*this)[Indices::canonicalToActiveComponentIndex(FluidSystem::waterCompIdx)] /=
175 FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, pvtRegionIdx);
177 if constexpr (enableSolvent) {
178 (*this)[Indices::contiSolventEqIdx] /=
179 solventPvt.referenceDensity(pvtRegionIdx);
187 template <
class Flu
idState,
class RhsEval>
190 const RhsEval& volume)
192 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
193 (*this)[conti0EqIdx + compIdx] =
194 fluidState.density(phaseIdx) *
195 fluidState.massFraction(phaseIdx, compIdx) *
203 template <
class RhsEval>
206 for (
unsigned i = 0; i < this->size(); ++i) {
Contains the classes required to extend the black-oil model by brine.
Contains the classes required to extend the black-oil model to include the effects of foam.
Contains the classes required to extend the black-oil model by MICP.
Contains the classes required to extend the black-oil model by polymer.
Contains the classes required to extend the black-oil model by solvents.
Contains the high level supplements required to extend the black oil model by brine.
Definition: blackoilbrinemodules.hh:55
Contains the high level supplements required to extend the black oil model to include the effects of ...
Definition: blackoilfoammodules.hh:58
Contains the high level supplements required to extend the black oil model by MICP.
Definition: blackoilmicpmodules.hh:54
Contains the high level supplements required to extend the black oil model by polymer.
Definition: blackoilpolymermodules.hh:64
Scalar molarMass() const
Definition: blackoilpolymermodules.hh:555
Implements a vector representing mass, molar or volumetric rates for the black oil model.
Definition: blackoilratevector.hh:62
BlackOilRateVector(Scalar value)
Definition: blackoilratevector.hh:95
void setMolarRate(const ParentType &value, unsigned pvtRegionIdx=0)
Set a molar rate of the conservation quantities.
Definition: blackoilratevector.hh:130
BlackOilRateVector & operator=(const RhsEval &value)
Assignment operator from a scalar or a function evaluation.
Definition: blackoilratevector.hh:204
void setVolumetricRate(const FluidState &fluidState, unsigned phaseIdx, const RhsEval &volume)
Set a volumetric rate of a phase.
Definition: blackoilratevector.hh:188
void setMassRate(const ParentType &value, unsigned pvtRegionIdx=0)
Set a mass rate of the conservation quantities.
Definition: blackoilratevector.hh:101
BlackOilRateVector()
Definition: blackoilratevector.hh:89
Contains the high level supplements required to extend the black oil model by solvents.
Definition: blackoilsolventmodules.hh:66
static const SolventPvt & solventPvt()
Definition: blackoilsolventmodules.hh:363
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