DeadOilPvt.hpp
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*/
27#ifndef OPM_DEAD_OIL_PVT_HPP
28#define OPM_DEAD_OIL_PVT_HPP
29
31
32#if HAVE_ECL_INPUT
33#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
34#include <opm/input/eclipse/EclipseState/Tables/PvdoTable.hpp>
35#include <opm/input/eclipse/EclipseState/Tables/TableManager.hpp>
36#endif
37
38namespace Opm {
43template <class Scalar>
45{
46public:
48
49 DeadOilPvt() = default;
50 DeadOilPvt(const std::vector<Scalar>& oilReferenceDensity,
51 const std::vector<TabulatedOneDFunction>& inverseOilB,
52 const std::vector<TabulatedOneDFunction>& oilMu,
53 const std::vector<TabulatedOneDFunction>& inverseOilBMu)
54 : oilReferenceDensity_(oilReferenceDensity)
55 , inverseOilB_(inverseOilB)
56 , oilMu_(oilMu)
57 , inverseOilBMu_(inverseOilBMu)
58 { }
59
60#if HAVE_ECL_INPUT
64 void initFromState(const EclipseState& eclState, const Schedule&)
65 {
66 const auto& pvdoTables = eclState.getTableManager().getPvdoTables();
67 const auto& densityTable = eclState.getTableManager().getDensityTable();
68
69 assert(pvdoTables.size() == densityTable.size());
70
71 size_t numRegions = pvdoTables.size();
73
74 for (unsigned regionIdx = 0; regionIdx < numRegions; ++ regionIdx) {
75 Scalar rhoRefO = densityTable[regionIdx].oil;
76 Scalar rhoRefG = densityTable[regionIdx].gas;
77 Scalar rhoRefW = densityTable[regionIdx].water;
78
79 setReferenceDensities(regionIdx, rhoRefO, rhoRefG, rhoRefW);
80
81 const auto& pvdoTable = pvdoTables.getTable<PvdoTable>(regionIdx);
82
83 const auto& BColumn(pvdoTable.getFormationFactorColumn());
84 std::vector<Scalar> invBColumn(BColumn.size());
85 for (unsigned i = 0; i < invBColumn.size(); ++i)
86 invBColumn[i] = 1/BColumn[i];
87
88 inverseOilB_[regionIdx].setXYArrays(pvdoTable.numRows(),
89 pvdoTable.getPressureColumn(),
90 invBColumn);
91 oilMu_[regionIdx].setXYArrays(pvdoTable.numRows(),
92 pvdoTable.getPressureColumn(),
93 pvdoTable.getViscosityColumn());
94 }
95
96 initEnd();
97 }
98#endif // HAVE_ECL_INPUT
99
101 {
102 oilReferenceDensity_.resize(numRegions);
103 inverseOilB_.resize(numRegions);
104 inverseOilBMu_.resize(numRegions);
105 oilMu_.resize(numRegions);
106 }
107
111 void setReferenceDensities(unsigned regionIdx,
112 Scalar rhoRefOil,
113 Scalar /*rhoRefGas*/,
114 Scalar /*rhoRefWater*/)
115 {
116 oilReferenceDensity_[regionIdx] = rhoRefOil;
117 }
118
129 void setInverseOilFormationVolumeFactor(unsigned regionIdx, const TabulatedOneDFunction& invBo)
130 { inverseOilB_[regionIdx] = invBo; }
131
137 void setOilViscosity(unsigned regionIdx, const TabulatedOneDFunction& muo)
138 { oilMu_[regionIdx] = muo; }
139
143 void initEnd()
144 {
145 // calculate the final 2D functions which are used for interpolation.
146 size_t numRegions = oilMu_.size();
147 for (unsigned regionIdx = 0; regionIdx < numRegions; ++ regionIdx) {
148 // calculate the table which stores the inverse of the product of the oil
149 // formation volume factor and the oil viscosity
150 const auto& oilMu = oilMu_[regionIdx];
151 const auto& invOilB = inverseOilB_[regionIdx];
152 assert(oilMu.numSamples() == invOilB.numSamples());
153
154 std::vector<Scalar> invBMuColumn;
155 std::vector<Scalar> pressureColumn;
156 invBMuColumn.resize(oilMu.numSamples());
157 pressureColumn.resize(oilMu.numSamples());
158
159 for (unsigned pIdx = 0; pIdx < oilMu.numSamples(); ++pIdx) {
160 pressureColumn[pIdx] = invOilB.xAt(pIdx);
161 invBMuColumn[pIdx] = invOilB.valueAt(pIdx)*1/oilMu.valueAt(pIdx);
162 }
163
164 inverseOilBMu_[regionIdx].setXYArrays(pressureColumn.size(),
165 pressureColumn,
166 invBMuColumn);
167 }
168 }
169
173 unsigned numRegions() const
174 { return inverseOilBMu_.size(); }
175
179 template <class Evaluation>
180 Evaluation internalEnergy(unsigned,
181 const Evaluation&,
182 const Evaluation&,
183 const Evaluation&) const
184 {
185 throw std::runtime_error("Requested the enthalpy of oil but the thermal option is not enabled");
186 }
187
191 template <class Evaluation>
192 Evaluation viscosity(unsigned regionIdx,
193 const Evaluation& temperature,
194 const Evaluation& pressure,
195 const Evaluation& /*Rs*/) const
196 { return saturatedViscosity(regionIdx, temperature, pressure); }
197
201 template <class Evaluation>
202 Evaluation saturatedViscosity(unsigned regionIdx,
203 const Evaluation& /*temperature*/,
204 const Evaluation& pressure) const
205 {
206 const Evaluation& invBo = inverseOilB_[regionIdx].eval(pressure, /*extrapolate=*/true);
207 const Evaluation& invMuoBo = inverseOilBMu_[regionIdx].eval(pressure, /*extrapolate=*/true);
208
209 return invBo/invMuoBo;
210 }
211
215 template <class Evaluation>
216 Evaluation inverseFormationVolumeFactor(unsigned regionIdx,
217 const Evaluation& /*temperature*/,
218 const Evaluation& pressure,
219 const Evaluation& /*Rs*/) const
220 { return inverseOilB_[regionIdx].eval(pressure, /*extrapolate=*/true); }
221
227 template <class Evaluation>
228 Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx,
229 const Evaluation& /*temperature*/,
230 const Evaluation& pressure) const
231 { return inverseOilB_[regionIdx].eval(pressure, /*extrapolate=*/true); }
232
236 template <class Evaluation>
237 Evaluation saturatedGasDissolutionFactor(unsigned /*regionIdx*/,
238 const Evaluation& /*temperature*/,
239 const Evaluation& /*pressure*/) const
240 { return 0.0; /* this is dead oil! */ }
241
245 template <class Evaluation>
246 Evaluation saturatedGasDissolutionFactor(unsigned /*regionIdx*/,
247 const Evaluation& /*temperature*/,
248 const Evaluation& /*pressure*/,
249 const Evaluation& /*oilSaturation*/,
250 const Evaluation& /*maxOilSaturation*/) const
251 { return 0.0; /* this is dead oil! */ }
252
259 template <class Evaluation>
260 Evaluation saturationPressure(unsigned /*regionIdx*/,
261 const Evaluation& /*temperature*/,
262 const Evaluation& /*Rs*/) const
263 { return 0.0; /* this is dead oil, so there isn't any meaningful saturation pressure! */ }
264
265 template <class Evaluation>
266 Evaluation saturatedGasMassFraction(unsigned /*regionIdx*/,
267 const Evaluation& /*temperature*/,
268 const Evaluation& /*pressure*/) const
269 { return 0.0; /* this is dead oil! */ }
270
271 template <class Evaluation>
272 Evaluation saturatedGasMoleFraction(unsigned /*regionIdx*/,
273 const Evaluation& /*temperature*/,
274 const Evaluation& /*pressure*/) const
275 { return 0.0; /* this is dead oil! */ }
276
277 template <class Evaluation>
278 Evaluation diffusionCoefficient(const Evaluation& /*temperature*/,
279 const Evaluation& /*pressure*/,
280 unsigned /*compIdx*/) const
281 {
282 throw std::runtime_error("Not implemented: The PVT model does not provide a diffusionCoefficient()");
283 }
284
285 const Scalar oilReferenceDensity(unsigned regionIdx) const
286 { return oilReferenceDensity_[regionIdx]; }
287
288 const std::vector<TabulatedOneDFunction>& inverseOilB() const
289 { return inverseOilB_; }
290
291 const std::vector<TabulatedOneDFunction>& oilMu() const
292 { return oilMu_; }
293
294 const std::vector<TabulatedOneDFunction>& inverseOilBMu() const
295 { return inverseOilBMu_; }
296
297 bool operator==(const DeadOilPvt<Scalar>& data) const
298 {
299 return this->oilReferenceDensity_ == data.oilReferenceDensity_ &&
300 this->inverseOilB() == data.inverseOilB() &&
301 this->oilMu() == data.oilMu() &&
302 this->inverseOilBMu() == data.inverseOilBMu();
303 }
304
305private:
306 std::vector<Scalar> oilReferenceDensity_;
307 std::vector<TabulatedOneDFunction> inverseOilB_;
308 std::vector<TabulatedOneDFunction> oilMu_;
309 std::vector<TabulatedOneDFunction> inverseOilBMu_;
310};
311
312} // namespace Opm
313
314#endif
This class represents the Pressure-Volume-Temperature relations of the oil phase without dissolved ga...
Definition: DeadOilPvt.hpp:45
bool operator==(const DeadOilPvt< Scalar > &data) const
Definition: DeadOilPvt.hpp:297
void setOilViscosity(unsigned regionIdx, const TabulatedOneDFunction &muo)
Initialize the viscosity of the oil phase.
Definition: DeadOilPvt.hpp:137
Evaluation saturatedGasDissolutionFactor(unsigned, const Evaluation &, const Evaluation &) const
Returns the gas dissolution factor [m^3/m^3] of the oil phase.
Definition: DeadOilPvt.hpp:237
void setInverseOilFormationVolumeFactor(unsigned regionIdx, const TabulatedOneDFunction &invBo)
Initialize the function for the oil formation volume factor.
Definition: DeadOilPvt.hpp:129
Evaluation viscosity(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &) const
Returns the dynamic viscosity [Pa s] of the fluid phase given a set of parameters.
Definition: DeadOilPvt.hpp:192
const std::vector< TabulatedOneDFunction > & inverseOilB() const
Definition: DeadOilPvt.hpp:288
Evaluation saturatedGasMassFraction(unsigned, const Evaluation &, const Evaluation &) const
Definition: DeadOilPvt.hpp:266
void setReferenceDensities(unsigned regionIdx, Scalar rhoRefOil, Scalar, Scalar)
Initialize the reference densities of all fluids for a given PVT region.
Definition: DeadOilPvt.hpp:111
void initEnd()
Finish initializing the oil phase PVT properties.
Definition: DeadOilPvt.hpp:143
DeadOilPvt()=default
Evaluation saturatedViscosity(unsigned regionIdx, const Evaluation &, const Evaluation &pressure) const
Returns the dynamic viscosity [Pa s] of gas saturated oil given a pressure.
Definition: DeadOilPvt.hpp:202
Evaluation diffusionCoefficient(const Evaluation &, const Evaluation &, unsigned) const
Definition: DeadOilPvt.hpp:278
const std::vector< TabulatedOneDFunction > & oilMu() const
Definition: DeadOilPvt.hpp:291
DeadOilPvt(const std::vector< Scalar > &oilReferenceDensity, const std::vector< TabulatedOneDFunction > &inverseOilB, const std::vector< TabulatedOneDFunction > &oilMu, const std::vector< TabulatedOneDFunction > &inverseOilBMu)
Definition: DeadOilPvt.hpp:50
const std::vector< TabulatedOneDFunction > & inverseOilBMu() const
Definition: DeadOilPvt.hpp:294
void setNumRegions(size_t numRegions)
Definition: DeadOilPvt.hpp:100
Evaluation saturatedGasMoleFraction(unsigned, const Evaluation &, const Evaluation &) const
Definition: DeadOilPvt.hpp:272
Evaluation internalEnergy(unsigned, const Evaluation &, const Evaluation &, const Evaluation &) const
Returns the specific enthalpy [J/kg] of oil given a set of parameters.
Definition: DeadOilPvt.hpp:180
Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &, const Evaluation &pressure) const
Returns the formation volume factor [-] of saturated oil.
Definition: DeadOilPvt.hpp:228
const Scalar oilReferenceDensity(unsigned regionIdx) const
Definition: DeadOilPvt.hpp:285
Evaluation saturationPressure(unsigned, const Evaluation &, const Evaluation &) const
Returns the saturation pressure of the oil phase [Pa] depending on its mass fraction of the gas compo...
Definition: DeadOilPvt.hpp:260
Evaluation inverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &, const Evaluation &pressure, const Evaluation &) const
Returns the formation volume factor [-] of the fluid phase.
Definition: DeadOilPvt.hpp:216
unsigned numRegions() const
Return the number of PVT regions which are considered by this PVT-object.
Definition: DeadOilPvt.hpp:173
Evaluation saturatedGasDissolutionFactor(unsigned, const Evaluation &, const Evaluation &, const Evaluation &, const Evaluation &) const
Returns the gas dissolution factor [m^3/m^3] of the oil phase.
Definition: DeadOilPvt.hpp:246
Implements a linearly interpolated scalar function that depends on one variable.
Definition: Tabulated1DFunction.hpp:48
Definition: Air_Mesitylene.hpp:34