blackoilratevector.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef EWOMS_BLACK_OIL_RATE_VECTOR_HH
29#define EWOMS_BLACK_OIL_RATE_VECTOR_HH
30
31#include <dune/common/fvector.hh>
32
33#include <opm/material/common/Valgrind.hpp>
34#include <opm/material/constraintsolvers/NcpFlash.hpp>
35
37
38namespace Opm {
39
49template <class TypeTag>
51 : public Dune::FieldVector<GetPropType<TypeTag, Properties::Evaluation>,
52 getPropValue<TypeTag, Properties::NumEq>()>
53{
58
64
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>;
78
79public:
80 BlackOilRateVector() : ParentType()
81 { Valgrind::SetUndefined(*this); }
82
86 BlackOilRateVector(Scalar value) : ParentType(Toolbox::createConstant(value))
87 {}
88
92 void setMassRate(const ParentType& value, unsigned pvtRegionIdx = 0)
93 {
94 ParentType::operator=(value);
95
96 // convert to "surface volume" if requested
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);
101 }
102 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
103 (*this)[Indices::canonicalToActiveComponentIndex(FluidSystem::oilCompIdx)] /=
104 FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, pvtRegionIdx);
105 }
106 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
107 (*this)[Indices::canonicalToActiveComponentIndex(FluidSystem::waterCompIdx)] /=
108 FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, pvtRegionIdx);
109 }
110 if constexpr (enableSolvent) {
111 const auto& solventPvt = SolventModule::solventPvt();
112 (*this)[Indices::contiSolventEqIdx] /=
113 solventPvt.referenceDensity(pvtRegionIdx);
114 }
115
116 }
117 }
118
122 void setMolarRate(const ParentType& value, unsigned pvtRegionIdx = 0)
123 {
124 // first, assign molar rates
125 ParentType::operator=(value);
126
127 // then, convert them to mass rates
128 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
129 (*this)[conti0EqIdx + compIdx] *= FluidSystem::molarMass(compIdx, pvtRegionIdx);
130
131 const auto& solventPvt = SolventModule::solventPvt();
132 (*this)[Indices::contiSolventEqIdx] *= solventPvt.molarMass(pvtRegionIdx);
133
134 if constexpr (enablePolymer) {
135 if constexpr (enablePolymerMolarWeight )
136 throw std::logic_error("Set molar rate with polymer weight tracking not implemented");
137
138 (*this)[Indices::contiPolymerEqIdx] *= PolymerModule::molarMass(pvtRegionIdx);
139 }
140
141 if constexpr (enableFoam) {
142 throw std::logic_error("setMolarRate() not implemented for foam");
143 }
144
145 if constexpr (enableBrine) {
146 throw std::logic_error("setMolarRate() not implemented for salt water");
147 }
148
149 if constexpr (enableMICP) {
150 throw std::logic_error("setMolarRate() not implemented for MICP");
151 }
152
153 // convert to "surface volume" if requested
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);
158 }
159 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
160 (*this)[Indices::canonicalToActiveComponentIndex(FluidSystem::oilCompIdx)] /=
161 FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, pvtRegionIdx);
162 }
163 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
164 (*this)[Indices::canonicalToActiveComponentIndex(FluidSystem::waterCompIdx)] /=
165 FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, pvtRegionIdx);
166 }
167 if constexpr (enableSolvent) {
168 (*this)[Indices::contiSolventEqIdx] /=
169 solventPvt.referenceDensity(pvtRegionIdx);
170 }
171 }
172 }
173
177 template <class FluidState, class RhsEval>
178 void setVolumetricRate(const FluidState& fluidState,
179 unsigned phaseIdx,
180 const RhsEval& volume)
181 {
182 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
183 (*this)[conti0EqIdx + compIdx] =
184 fluidState.density(phaseIdx)
185 * fluidState.massFraction(phaseIdx, compIdx)
186 * volume;
187 }
188
192 template <class RhsEval>
193 BlackOilRateVector& operator=(const RhsEval& value)
194 {
195 for (unsigned i=0; i < this->size(); ++i)
196 (*this)[i] = value;
197 return *this;
198 }
199};
200
201} // namespace Opm
202
203#endif
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:57
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:68
static const SolventPvt & solventPvt()
Definition: blackoilsolventmodules.hh:607
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:242