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
41
42#include <stdexcept>
43
44namespace Opm {
45
55template <class TypeTag>
57 : public Dune::FieldVector<GetPropType<TypeTag, Properties::Evaluation>,
58 getPropValue<TypeTag, Properties::NumEq>()>
59{
64
65 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
66 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
67 enum { conti0EqIdx = Indices::conti0EqIdx };
68
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>();
75
76 using Toolbox = MathToolbox<Evaluation>;
77 using ParentType = Dune::FieldVector<Evaluation, numEq>;
78
79 using BrineModule = BlackOilBrineModule<TypeTag, enableBrine>;
81 using PolymerModule = BlackOilPolymerModule<TypeTag, enablePolymer>;
82 using SolventModule = BlackOilSolventModule<TypeTag, enableSolvent>;
83
84public:
85 BlackOilRateVector() : ParentType()
86 { Valgrind::SetUndefined(*this); }
87
91 BlackOilRateVector(Scalar value) : ParentType(Toolbox::createConstant(value))
92 {}
93
100 void setMassRate(const ParentType& value, unsigned pvtRegionIdx = 0)
101 {
102 ParentType::operator=(value);
103
104 // convert to "surface volume" if requested
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);
109 }
110 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
111 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx)] /=
112 FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, pvtRegionIdx);
113 }
114 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
115 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx)] /=
116 FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, pvtRegionIdx);
117 }
118 if constexpr (enableSolvent) {
119 const auto& solventPvt = SolventModule::solventPvt();
120 (*this)[Indices::contiSolventEqIdx] /=
121 solventPvt.referenceDensity(pvtRegionIdx);
122 }
123 }
124 }
125
132 void setMolarRate(const ParentType& value, unsigned pvtRegionIdx = 0)
133 {
134 // first, assign molar rates
135 ParentType::operator=(value);
136
137 // then, convert them to mass rates
138 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
139 (*this)[conti0EqIdx + compIdx] *= FluidSystem::molarMass(compIdx, pvtRegionIdx);
140 }
141
142 if constexpr (enableSolvent) {
143 const auto& solventPvt = SolventModule::solventPvt();
144 (*this)[Indices::contiSolventEqIdx] *= solventPvt.molarMass(pvtRegionIdx);
145 }
146
147 if constexpr (enablePolymer) {
148 if constexpr (enablePolymerMolarWeight) {
149 throw std::logic_error("Set molar rate with polymer weight tracking not implemented");
150 }
151
152 (*this)[Indices::contiPolymerEqIdx] *= PolymerModule::molarMass(pvtRegionIdx);
153 }
154
155 if constexpr (enableFoam) {
156 throw std::logic_error("setMolarRate() not implemented for foam");
157 }
158
159 if constexpr (enableBrine) {
160 throw std::logic_error("setMolarRate() not implemented for salt water");
161 }
162
163 if constexpr (enableBioeffects) {
164 throw std::logic_error("setMolarRate() not implemented for bioeffects (biofilm/MICP)");
165 }
166
167 // convert to "surface volume" if requested
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);
172 }
173 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
174 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx)] /=
175 FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, pvtRegionIdx);
176 }
177 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
178 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx)] /=
179 FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, pvtRegionIdx);
180 }
181 if constexpr (enableSolvent) {
182 const auto& solventPvt = SolventModule::solventPvt();
183 (*this)[Indices::contiSolventEqIdx] /=
184 solventPvt.referenceDensity(pvtRegionIdx);
185 }
186 }
187 }
188
192 template <class FluidState, class RhsEval>
193 void setVolumetricRate(const FluidState& fluidState,
194 unsigned phaseIdx,
195 const RhsEval& volume)
196 {
197 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
198 (*this)[conti0EqIdx + compIdx] =
199 fluidState.density(phaseIdx) *
200 fluidState.massFraction(phaseIdx, compIdx) *
201 volume;
202 }
203 }
204
208 template <class RhsEval>
209 BlackOilRateVector& operator=(const RhsEval& value)
210 {
211 for (unsigned i = 0; i < this->size(); ++i) {
212 (*this)[i] = value;
213 }
214 return *this;
215 }
216};
217
218} // namespace Opm
219
220#endif
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