27 #ifndef OPM_DEAD_OIL_PVT_HPP
28 #define OPM_DEAD_OIL_PVT_HPP
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>
43 template <
class Scalar>
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)
57 , inverseOilBMu_(inverseOilBMu)
64 void initFromState(
const EclipseState& eclState,
const Schedule&)
66 const auto& pvdoTables = eclState.getTableManager().getPvdoTables();
67 const auto& densityTable = eclState.getTableManager().getDensityTable();
69 assert(pvdoTables.size() == densityTable.size());
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;
81 const auto& pvdoTable = pvdoTables.getTable<PvdoTable>(regionIdx);
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];
88 inverseOilB_[regionIdx].setXYArrays(pvdoTable.numRows(),
89 pvdoTable.getPressureColumn(),
91 oilMu_[regionIdx].setXYArrays(pvdoTable.numRows(),
92 pvdoTable.getPressureColumn(),
93 pvdoTable.getViscosityColumn());
116 oilReferenceDensity_[regionIdx] = rhoRefOil;
130 { inverseOilB_[regionIdx] = invBo; }
138 { oilMu_[regionIdx] = muo; }
147 for (
unsigned regionIdx = 0; regionIdx <
numRegions; ++ regionIdx) {
150 const auto& oilMu = oilMu_[regionIdx];
151 const auto& invOilB = inverseOilB_[regionIdx];
152 assert(oilMu.numSamples() == invOilB.numSamples());
154 std::vector<Scalar> invBMuColumn;
155 std::vector<Scalar> pressureColumn;
156 invBMuColumn.resize(oilMu.numSamples());
157 pressureColumn.resize(oilMu.numSamples());
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);
164 inverseOilBMu_[regionIdx].setXYArrays(pressureColumn.size(),
174 {
return inverseOilBMu_.size(); }
179 template <
class Evaluation>
183 const Evaluation&)
const
185 throw std::runtime_error(
"Requested the enthalpy of oil but the thermal option is not enabled");
191 template <
class Evaluation>
193 const Evaluation& temperature,
194 const Evaluation& pressure,
195 const Evaluation& )
const
201 template <
class Evaluation>
204 const Evaluation& pressure)
const
206 const Evaluation& invBo = inverseOilB_[regionIdx].eval(pressure,
true);
207 const Evaluation& invMuoBo = inverseOilBMu_[regionIdx].eval(pressure,
true);
209 return invBo/invMuoBo;
215 template <
class Evaluation>
218 const Evaluation& pressure,
219 const Evaluation& )
const
220 {
return inverseOilB_[regionIdx].eval(pressure,
true); }
227 template <
class Evaluation>
230 const Evaluation& pressure)
const
231 {
return inverseOilB_[regionIdx].eval(pressure,
true); }
236 template <
class Evaluation>
239 const Evaluation& )
const
245 template <
class Evaluation>
250 const Evaluation& )
const
259 template <
class Evaluation>
262 const Evaluation& )
const
265 template <
class Evaluation>
266 Evaluation saturatedGasMassFraction(
unsigned ,
268 const Evaluation& )
const
271 template <
class Evaluation>
272 Evaluation saturatedGasMoleFraction(
unsigned ,
274 const Evaluation& )
const
277 template <
class Evaluation>
278 Evaluation diffusionCoefficient(
const Evaluation& ,
282 throw std::runtime_error(
"Not implemented: The PVT model does not provide a diffusionCoefficient()");
285 const Scalar oilReferenceDensity(
unsigned regionIdx)
const
286 {
return oilReferenceDensity_[regionIdx]; }
288 const std::vector<TabulatedOneDFunction>& inverseOilB()
const
289 {
return inverseOilB_; }
291 const std::vector<TabulatedOneDFunction>& oilMu()
const
294 const std::vector<TabulatedOneDFunction>& inverseOilBMu()
const
295 {
return inverseOilBMu_; }
297 bool operator==(
const DeadOilPvt<Scalar>& data)
const
299 return this->oilReferenceDensity_ == data.oilReferenceDensity_ &&
300 this->inverseOilB() == data.inverseOilB() &&
301 this->oilMu() == data.oilMu() &&
302 this->inverseOilBMu() == data.inverseOilBMu();
306 std::vector<Scalar> oilReferenceDensity_;
307 std::vector<TabulatedOneDFunction> inverseOilB_;
308 std::vector<TabulatedOneDFunction> oilMu_;
309 std::vector<TabulatedOneDFunction> inverseOilBMu_;
Implements a linearly interpolated scalar function that depends on one variable.
This class represents the Pressure-Volume-Temperature relations of the oil phase without dissolved ga...
Definition: DeadOilPvt.hpp:45
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
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
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 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
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