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>
59class BlackOilRateVector
60 :
public Dune::FieldVector<GetPropType<TypeTag, Properties::Evaluation>,
61 getPropValue<TypeTag, Properties::NumEq>()>
75 enum { conti0EqIdx = Indices::conti0EqIdx };
82 using Toolbox = MathToolbox<Evaluation>;
83 using ParentType = Dune::FieldVector<Evaluation, numEq>;
86 BlackOilRateVector() : ParentType()
87 { Valgrind::SetUndefined(*
this); }
101 void setMassRate(
const ParentType& value,
unsigned pvtRegionIdx = 0)
103 ParentType::operator=(value);
107 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
108 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx)] /=
109 FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, pvtRegionIdx);
111 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
112 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx)] /=
113 FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, pvtRegionIdx);
115 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
116 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx)] /=
117 FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, pvtRegionIdx);
119 if constexpr (enableSolvent) {
120 const auto& solventPvt = SolventModule::solventPvt();
121 (*this)[Indices::contiSolventEqIdx] /=
122 solventPvt.referenceDensity(pvtRegionIdx);
136 ParentType::operator=(value);
139 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
140 (*this)[conti0EqIdx + compIdx] *= FluidSystem::molarMass(compIdx, pvtRegionIdx);
143 const auto& solventPvt = SolventModule::solventPvt();
144 (*this)[Indices::contiSolventEqIdx] *= solventPvt.molarMass(pvtRegionIdx);
146 if constexpr (enablePolymer) {
147 if constexpr (enablePolymerMolarWeight) {
148 throw std::logic_error(
"Set molar rate with polymer weight tracking not implemented");
151 (*this)[Indices::contiPolymerEqIdx] *= PolymerModule::molarMass(pvtRegionIdx);
154 if constexpr (enableFoam) {
155 throw std::logic_error(
"setMolarRate() not implemented for foam");
158 if constexpr (enableBrine) {
159 throw std::logic_error(
"setMolarRate() not implemented for salt water");
162 if constexpr (enableBioeffects) {
163 throw std::logic_error(
"setMolarRate() not implemented for bioeffects (biofilm/MICP)");
168 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
169 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx)] /=
170 FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, pvtRegionIdx);
172 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
173 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx)] /=
174 FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, pvtRegionIdx);
176 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
177 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx)] /=
178 FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, pvtRegionIdx);
180 if constexpr (enableSolvent) {
181 (*this)[Indices::contiSolventEqIdx] /=
182 solventPvt.referenceDensity(pvtRegionIdx);
190 template <
class Flu
idState,
class RhsEval>
193 const RhsEval& volume)
195 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
196 (*this)[conti0EqIdx + compIdx] =
197 fluidState.density(phaseIdx) *
198 fluidState.massFraction(phaseIdx, compIdx) *
206 template <
class RhsEval>
209 for (
unsigned i = 0; i < this->size(); ++i) {
Contains the classes required to extend the black-oil model by bioeffects.
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 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:58
Contains the high level supplements required to extend the black oil model to include the effects of ...
Definition blackoilfoammodules.hh:60
Contains the high level supplements required to extend the black oil model by polymer.
Definition blackoilpolymermodules.hh:65
BlackOilRateVector(Scalar value)
Definition blackoilratevector.hh:92
void setMolarRate(const ParentType &value, unsigned pvtRegionIdx=0)
Set a molar rate of the conservation quantities.
Definition blackoilratevector.hh:133
BlackOilRateVector & operator=(const RhsEval &value)
Assignment operator from a scalar or a function evaluation.
Definition blackoilratevector.hh:207
void setVolumetricRate(const FluidState &fluidState, unsigned phaseIdx, const RhsEval &volume)
Set a volumetric rate of a phase.
Definition blackoilratevector.hh:191
void setMassRate(const ParentType &value, unsigned pvtRegionIdx=0)
Set a mass rate of the conservation quantities.
Definition blackoilratevector.hh:101
Contains the high level supplements required to extend the black oil model by solvents.
Definition blackoilsolventmodules.hh:69
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.
This file contains a set of helper functions used by VFPProd / VFPInj.
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
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240