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>
55template <
class TypeTag>
57 :
public Dune::FieldVector<GetPropType<TypeTag, Properties::Evaluation>,
58 getPropValue<TypeTag, Properties::NumEq>()>
65 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
66 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
67 enum { conti0EqIdx = Indices::conti0EqIdx };
69 static constexpr bool enableBioeffects = getPropValue<TypeTag, Properties::EnableBioeffects>();
70 static constexpr bool enableBrine = getPropValue<TypeTag, Properties::EnableBrine>();
71 static constexpr bool enableFoam = getPropValue<TypeTag, Properties::EnableFoam>();
72 static constexpr bool enablePolymer = getPropValue<TypeTag, Properties::EnablePolymer>();
73 static constexpr bool enablePolymerMolarWeight = getPropValue<TypeTag, Properties::EnablePolymerMW>();
74 static constexpr bool enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>();
76 using Toolbox = MathToolbox<Evaluation>;
77 using ParentType = Dune::FieldVector<Evaluation, numEq>;
79 using BrineModule = BlackOilBrineModule<TypeTag, enableBrine>;
81 using PolymerModule = BlackOilPolymerModule<TypeTag, enablePolymer>;
82 using SolventModule = BlackOilSolventModule<TypeTag, enableSolvent>;
86 { Valgrind::SetUndefined(*
this); }
100 void setMassRate(
const ParentType& value,
unsigned pvtRegionIdx = 0)
102 ParentType::operator=(value);
105 if constexpr (getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>()) {
106 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
107 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx)] /=
108 FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, pvtRegionIdx);
110 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
111 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx)] /=
112 FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, pvtRegionIdx);
114 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
115 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx)] /=
116 FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, pvtRegionIdx);
118 if constexpr (enableSolvent) {
119 const auto& solventPvt = SolventModule::solventPvt();
120 (*this)[Indices::contiSolventEqIdx] /=
121 solventPvt.referenceDensity(pvtRegionIdx);
135 ParentType::operator=(value);
138 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
139 (*this)[conti0EqIdx + compIdx] *= FluidSystem::molarMass(compIdx, pvtRegionIdx);
142 if constexpr (enableSolvent) {
143 const auto& solventPvt = SolventModule::solventPvt();
144 (*this)[Indices::contiSolventEqIdx] *= solventPvt.molarMass(pvtRegionIdx);
147 if constexpr (enablePolymer) {
148 if constexpr (enablePolymerMolarWeight) {
149 throw std::logic_error(
"Set molar rate with polymer weight tracking not implemented");
152 (*this)[Indices::contiPolymerEqIdx] *= PolymerModule::molarMass(pvtRegionIdx);
155 if constexpr (enableFoam) {
156 throw std::logic_error(
"setMolarRate() not implemented for foam");
159 if constexpr (enableBrine) {
160 throw std::logic_error(
"setMolarRate() not implemented for salt water");
163 if constexpr (enableBioeffects) {
164 throw std::logic_error(
"setMolarRate() not implemented for bioeffects (biofilm/MICP)");
168 if constexpr (getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>()) {
169 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
170 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx)] /=
171 FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, pvtRegionIdx);
173 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
174 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx)] /=
175 FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, pvtRegionIdx);
177 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
178 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx)] /=
179 FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, pvtRegionIdx);
181 if constexpr (enableSolvent) {
182 const auto& solventPvt = SolventModule::solventPvt();
183 (*this)[Indices::contiSolventEqIdx] /=
184 solventPvt.referenceDensity(pvtRegionIdx);
192 template <
class Flu
idState,
class RhsEval>
195 const RhsEval& volume)
197 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
198 (*this)[conti0EqIdx + compIdx] =
199 fluidState.density(phaseIdx) *
200 fluidState.massFraction(phaseIdx, compIdx) *
208 template <
class RhsEval>
211 for (
unsigned i = 0; i < this->size(); ++i) {
Contains classes extending the black-oil model. \detail This file holds dummy definitions,...
Declares the properties required by the black oil model.
Definition: blackoilfoammodules.hh:57
Implements a vector representing mass, molar or volumetric rates for the black oil model.
Definition: blackoilratevector.hh:59
BlackOilRateVector(Scalar value)
Definition: blackoilratevector.hh:91
void setMolarRate(const ParentType &value, unsigned pvtRegionIdx=0)
Set a molar rate of the conservation quantities.
Definition: blackoilratevector.hh:132
BlackOilRateVector & operator=(const RhsEval &value)
Assignment operator from a scalar or a function evaluation.
Definition: blackoilratevector.hh:209
void setVolumetricRate(const FluidState &fluidState, unsigned phaseIdx, const RhsEval &volume)
Set a volumetric rate of a phase.
Definition: blackoilratevector.hh:193
void setMassRate(const ParentType &value, unsigned pvtRegionIdx=0)
Set a mass rate of the conservation quantities.
Definition: blackoilratevector.hh:100
BlackOilRateVector()
Definition: blackoilratevector.hh:85
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: blackoilbioeffectsmodules.hh:45
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