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
72
73 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
74 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
75 enum { conti0EqIdx = Indices::conti0EqIdx };
76 enum { enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>() };
77 enum { enablePolymer = getPropValue<TypeTag, Properties::EnablePolymer>() };
78 enum { enablePolymerMolarWeight = getPropValue<TypeTag, Properties::EnablePolymerMW>() };
79 enum { enableFoam = getPropValue<TypeTag, Properties::EnableFoam>() };
80 enum { enableBrine = getPropValue<TypeTag, Properties::EnableBrine>() };
81 enum { enableBioeffects = getPropValue<TypeTag, Properties::EnableBioeffects>() };
82 using Toolbox = MathToolbox<Evaluation>;
83 using ParentType = Dune::FieldVector<Evaluation, numEq>;
84
85public:
86 BlackOilRateVector() : ParentType()
87 { Valgrind::SetUndefined(*this); }
88
92 BlackOilRateVector(Scalar value) : ParentType(Toolbox::createConstant(value))
93 {}
94
101 void setMassRate(const ParentType& value, unsigned pvtRegionIdx = 0)
102 {
103 ParentType::operator=(value);
104
105 // convert to "surface volume" if requested
106 if constexpr (getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>()) {
107 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
108 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx)] /=
109 FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, pvtRegionIdx);
110 }
111 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
112 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx)] /=
113 FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, pvtRegionIdx);
114 }
115 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
116 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx)] /=
117 FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, pvtRegionIdx);
118 }
119 if constexpr (enableSolvent) {
120 const auto& solventPvt = SolventModule::solventPvt();
121 (*this)[Indices::contiSolventEqIdx] /=
122 solventPvt.referenceDensity(pvtRegionIdx);
123 }
124 }
125 }
126
133 void setMolarRate(const ParentType& value, unsigned pvtRegionIdx = 0)
134 {
135 // first, assign molar rates
136 ParentType::operator=(value);
137
138 // then, convert them to mass rates
139 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
140 (*this)[conti0EqIdx + compIdx] *= FluidSystem::molarMass(compIdx, pvtRegionIdx);
141 }
142
143 const auto& solventPvt = SolventModule::solventPvt();
144 (*this)[Indices::contiSolventEqIdx] *= solventPvt.molarMass(pvtRegionIdx);
145
146 if constexpr (enablePolymer) {
147 if constexpr (enablePolymerMolarWeight) {
148 throw std::logic_error("Set molar rate with polymer weight tracking not implemented");
149 }
150
151 (*this)[Indices::contiPolymerEqIdx] *= PolymerModule::molarMass(pvtRegionIdx);
152 }
153
154 if constexpr (enableFoam) {
155 throw std::logic_error("setMolarRate() not implemented for foam");
156 }
157
158 if constexpr (enableBrine) {
159 throw std::logic_error("setMolarRate() not implemented for salt water");
160 }
161
162 if constexpr (enableBioeffects) {
163 throw std::logic_error("setMolarRate() not implemented for bioeffects (biofilm/MICP)");
164 }
165
166 // convert to "surface volume" if requested
167 if constexpr (getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>()) {
168 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
169 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx)] /=
170 FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, pvtRegionIdx);
171 }
172 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
173 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx)] /=
174 FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, pvtRegionIdx);
175 }
176 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
177 (*this)[FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx)] /=
178 FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, pvtRegionIdx);
179 }
180 if constexpr (enableSolvent) {
181 (*this)[Indices::contiSolventEqIdx] /=
182 solventPvt.referenceDensity(pvtRegionIdx);
183 }
184 }
185 }
186
190 template <class FluidState, class RhsEval>
191 void setVolumetricRate(const FluidState& fluidState,
192 unsigned phaseIdx,
193 const RhsEval& volume)
194 {
195 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
196 (*this)[conti0EqIdx + compIdx] =
197 fluidState.density(phaseIdx) *
198 fluidState.massFraction(phaseIdx, compIdx) *
199 volume;
200 }
201 }
202
206 template <class RhsEval>
207 BlackOilRateVector& operator=(const RhsEval& value)
208 {
209 for (unsigned i = 0; i < this->size(); ++i) {
210 (*this)[i] = value;
211 }
212 return *this;
213 }
214};
215
216} // namespace Opm
217
218#endif
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: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 polymer.
Definition: blackoilpolymermodules.hh:64
Scalar molarMass() const
Definition: blackoilpolymermodules.hh:553
Implements a vector representing mass, molar or volumetric rates for the black oil model.
Definition: blackoilratevector.hh:62
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
BlackOilRateVector()
Definition: blackoilratevector.hh:86
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:363
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:43
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