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/MathToolbox.hpp>
34#include <opm/material/common/Valgrind.hpp>
35#include <opm/material/constraintsolvers/NcpFlash.hpp>
36
44
45#include <stdexcept>
46
47namespace Opm {
48
58template <class TypeTag>
60 : public Dune::FieldVector<GetPropType<TypeTag, Properties::Evaluation>,
61 getPropValue<TypeTag, Properties::NumEq>()>
62{
67
73
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>;
87
88public:
89 BlackOilRateVector() : ParentType()
90 { Valgrind::SetUndefined(*this); }
91
95 BlackOilRateVector(Scalar value) : ParentType(Toolbox::createConstant(value))
96 {}
97
104 void setMassRate(const ParentType& value, unsigned pvtRegionIdx = 0)
105 {
106 ParentType::operator=(value);
107
108 // convert to "surface volume" if requested
109 if constexpr (getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>()) {
110 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
111 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx)] /=
112 FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, pvtRegionIdx);
113 }
114 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
115 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx)] /=
116 FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, pvtRegionIdx);
117 }
118 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
119 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx)] /=
120 FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, pvtRegionIdx);
121 }
122 if constexpr (enableSolvent) {
123 const auto& solventPvt = SolventModule::solventPvt();
124 (*this)[Indices::contiSolventEqIdx] /=
125 solventPvt.referenceDensity(pvtRegionIdx);
126 }
127 }
128 }
129
136 void setMolarRate(const ParentType& value, unsigned pvtRegionIdx = 0)
137 {
138 // first, assign molar rates
139 ParentType::operator=(value);
140
141 // then, convert them to mass rates
142 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
143 (*this)[conti0EqIdx + compIdx] *= FluidSystem::molarMass(compIdx, pvtRegionIdx);
144 }
145
146 const auto& solventPvt = SolventModule::solventPvt();
147 (*this)[Indices::contiSolventEqIdx] *= solventPvt.molarMass(pvtRegionIdx);
148
149 if constexpr (enablePolymer) {
150 if constexpr (enablePolymerMolarWeight) {
151 throw std::logic_error("Set molar rate with polymer weight tracking not implemented");
152 }
153
154 (*this)[Indices::contiPolymerEqIdx] *= PolymerModule::molarMass(pvtRegionIdx);
155 }
156
157 if constexpr (enableFoam) {
158 throw std::logic_error("setMolarRate() not implemented for foam");
159 }
160
161 if constexpr (enableBrine) {
162 throw std::logic_error("setMolarRate() not implemented for salt water");
163 }
164
165 if constexpr (enableMICP) {
166 throw std::logic_error("setMolarRate() not implemented for MICP");
167 }
168
169 // convert to "surface volume" if requested
170 if constexpr (getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>()) {
171 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
172 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx)] /=
173 FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, pvtRegionIdx);
174 }
175 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
176 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx)] /=
177 FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, pvtRegionIdx);
178 }
179 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
180 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx)] /=
181 FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, pvtRegionIdx);
182 }
183 if constexpr (enableSolvent) {
184 (*this)[Indices::contiSolventEqIdx] /=
185 solventPvt.referenceDensity(pvtRegionIdx);
186 }
187 }
188 }
189
193 template <class FluidState, class RhsEval>
194 void setVolumetricRate(const FluidState& fluidState,
195 unsigned phaseIdx,
196 const RhsEval& volume)
197 {
198 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
199 (*this)[conti0EqIdx + compIdx] =
200 fluidState.density(phaseIdx) *
201 fluidState.massFraction(phaseIdx, compIdx) *
202 volume;
203 }
204 }
205
209 template <class RhsEval>
210 BlackOilRateVector& operator=(const RhsEval& value)
211 {
212 for (unsigned i = 0; i < this->size(); ++i) {
213 (*this)[i] = value;
214 }
215 return *this;
216 }
217};
218
219} // namespace Opm
220
221#endif
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:56
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:136
BlackOilRateVector & operator=(const RhsEval &value)
Assignment operator from a scalar or a function evaluation.
Definition: blackoilratevector.hh:210
void setVolumetricRate(const FluidState &fluidState, unsigned phaseIdx, const RhsEval &volume)
Set a volumetric rate of a phase.
Definition: blackoilratevector.hh:194
void setMassRate(const ParentType &value, unsigned pvtRegionIdx=0)
Set a mass rate of the conservation quantities.
Definition: blackoilratevector.hh:104
BlackOilRateVector()
Definition: blackoilratevector.hh:89
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:365
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