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/Valgrind.hpp>
34#include <opm/material/constraintsolvers/NcpFlash.hpp>
49template <
class TypeTag>
51 :
public Dune::FieldVector<GetPropType<TypeTag, Properties::Evaluation>,
52 getPropValue<TypeTag, Properties::NumEq>()>
65 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
66 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
67 enum { conti0EqIdx = Indices::conti0EqIdx };
68 enum { contiEnergyEqIdx = Indices::contiEnergyEqIdx };
69 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
70 enum { enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>() };
71 enum { enablePolymer = getPropValue<TypeTag, Properties::EnablePolymer>() };
72 enum { enablePolymerMolarWeight = getPropValue<TypeTag, Properties::EnablePolymerMW>() };
73 enum { enableFoam = getPropValue<TypeTag, Properties::EnableFoam>() };
74 enum { enableBrine = getPropValue<TypeTag, Properties::EnableBrine>() };
75 enum { enableMICP = getPropValue<TypeTag, Properties::EnableMICP>() };
76 using Toolbox = MathToolbox<Evaluation>;
77 using ParentType = Dune::FieldVector<Evaluation, numEq>;
81 { Valgrind::SetUndefined(*
this); }
92 void setMassRate(
const ParentType& value,
unsigned pvtRegionIdx = 0)
94 ParentType::operator=(value);
97 if constexpr (getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>()) {
98 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
99 (*this)[Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx)] /=
100 FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, pvtRegionIdx);
102 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
103 (*this)[Indices::canonicalToActiveComponentIndex(FluidSystem::oilCompIdx)] /=
104 FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, pvtRegionIdx);
106 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
107 (*this)[Indices::canonicalToActiveComponentIndex(FluidSystem::waterCompIdx)] /=
108 FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, pvtRegionIdx);
110 if constexpr (enableSolvent) {
112 (*this)[Indices::contiSolventEqIdx] /=
113 solventPvt.referenceDensity(pvtRegionIdx);
125 ParentType::operator=(value);
128 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
129 (*
this)[conti0EqIdx + compIdx] *= FluidSystem::molarMass(compIdx, pvtRegionIdx);
132 (*this)[Indices::contiSolventEqIdx] *= solventPvt.molarMass(pvtRegionIdx);
134 if constexpr (enablePolymer) {
135 if constexpr (enablePolymerMolarWeight )
136 throw std::logic_error(
"Set molar rate with polymer weight tracking not implemented");
141 if constexpr (enableFoam) {
142 throw std::logic_error(
"setMolarRate() not implemented for foam");
145 if constexpr (enableBrine) {
146 throw std::logic_error(
"setMolarRate() not implemented for salt water");
149 if constexpr (enableMICP) {
150 throw std::logic_error(
"setMolarRate() not implemented for MICP");
154 if constexpr (getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>()) {
155 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
156 (*this)[Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx)] /=
157 FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, pvtRegionIdx);
159 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
160 (*this)[Indices::canonicalToActiveComponentIndex(FluidSystem::oilCompIdx)] /=
161 FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, pvtRegionIdx);
163 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
164 (*this)[Indices::canonicalToActiveComponentIndex(FluidSystem::waterCompIdx)] /=
165 FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, pvtRegionIdx);
167 if constexpr (enableSolvent) {
168 (*this)[Indices::contiSolventEqIdx] /=
169 solventPvt.referenceDensity(pvtRegionIdx);
177 template <
class Flu
idState,
class RhsEval>
180 const RhsEval& volume)
182 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
183 (*
this)[conti0EqIdx + compIdx] =
184 fluidState.density(phaseIdx)
185 * fluidState.massFraction(phaseIdx, compIdx)
192 template <
class RhsEval>
195 for (
unsigned i=0; i < this->size(); ++i)
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:59
Contains the high level supplements required to extend the black oil model by MICP.
Definition: blackoilmicpmodules.hh:56
Contains the high level supplements required to extend the black oil model by polymer.
Definition: blackoilpolymermodules.hh:61
const Scalar molarMass() const
Definition: blackoilpolymermodules.hh:784
Implements a vector representing mass, molar or volumetric rates for the black oil model.
Definition: blackoilratevector.hh:53
BlackOilRateVector(Scalar value)
Definition: blackoilratevector.hh:86
void setMolarRate(const ParentType &value, unsigned pvtRegionIdx=0)
Set a molar rate of the conservation quantities.
Definition: blackoilratevector.hh:122
BlackOilRateVector & operator=(const RhsEval &value)
Assignment operator from a scalar or a function evaluation.
Definition: blackoilratevector.hh:193
void setVolumetricRate(const FluidState &fluidState, unsigned phaseIdx, const RhsEval &volume)
Set a volumetric rate of a phase.
Definition: blackoilratevector.hh:178
void setMassRate(const ParentType &value, unsigned pvtRegionIdx=0)
Set a mass rate of the conservation quantities.
Definition: blackoilratevector.hh:92
BlackOilRateVector()
Definition: blackoilratevector.hh:80
Contains the high level supplements required to extend the black oil model by solvents.
Definition: blackoilsolventmodules.hh:69
static const SolventPvt & solventPvt()
Definition: blackoilsolventmodules.hh:608
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