mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
Add support for PBDV and PDDV in initstateequil.hh
This commit is contained in:
@@ -241,6 +241,138 @@ private:
|
||||
}
|
||||
};
|
||||
|
||||
/**
|
||||
* Type that implements "dissolved gas-oil ratio"
|
||||
* tabulated as a function of depth policy. Data
|
||||
* typically from keyword 'PBVD'.
|
||||
*/
|
||||
template <class FluidSystem>
|
||||
class PBVD : public RsFunction
|
||||
{
|
||||
public:
|
||||
/**
|
||||
* Constructor.
|
||||
*
|
||||
* \param[in] pvtRegionIdx The pvt region index
|
||||
* \param[in] depth Depth nodes.
|
||||
* \param[in] pbub Bubble-point pressure at @c depth.
|
||||
*/
|
||||
PBVD(const int pvtRegionIdx,
|
||||
const std::vector<double>& depth,
|
||||
const std::vector<double>& pbub)
|
||||
: pvtRegionIdx_(pvtRegionIdx)
|
||||
, pbubVsDepth_(depth, pbub)
|
||||
{}
|
||||
|
||||
/**
|
||||
* Function call.
|
||||
*
|
||||
* \param[in] depth Depth at which to calculate RS
|
||||
* value.
|
||||
*
|
||||
* \param[in] Pressure in the cell
|
||||
*
|
||||
* \param[in] temp Temperature at which to calculate RS
|
||||
* value.
|
||||
*
|
||||
* \return Dissolved gas-oil ratio (RS) at depth @c
|
||||
* depth and pressure @c press.
|
||||
*/
|
||||
double operator()(const double depth,
|
||||
const double cellPress,
|
||||
const double temp,
|
||||
const double satGas = 0.0) const
|
||||
{
|
||||
double press = cellPress;
|
||||
if (satGas <= 0.0) {
|
||||
if (pbubVsDepth_.xMin() > depth)
|
||||
press = pbubVsDepth_.valueAt(0);
|
||||
else if (pbubVsDepth_.xMax() < depth)
|
||||
press = pbubVsDepth_.valueAt(pbubVsDepth_.numSamples() - 1);
|
||||
else
|
||||
press = pbubVsDepth_.eval(depth, /*extrapolate=*/false);
|
||||
}
|
||||
return satRs(std::min(press, cellPress), temp);
|
||||
}
|
||||
|
||||
private:
|
||||
typedef Opm::Tabulated1DFunction<double> PbubVsDepthFunc;
|
||||
|
||||
const int pvtRegionIdx_;
|
||||
PbubVsDepthFunc pbubVsDepth_;
|
||||
|
||||
double satRs(const double press, const double temp) const
|
||||
{
|
||||
return FluidSystem::oilPvt().saturatedGasDissolutionFactor(pvtRegionIdx_, temp, press);
|
||||
}
|
||||
};
|
||||
|
||||
/**
|
||||
* Type that implements "vaporized oil-gas ratio"
|
||||
* tabulated as a function of depth policy. Data
|
||||
* taken from keyword 'PDVD'.
|
||||
*/
|
||||
template <class FluidSystem>
|
||||
class PDVD : public RsFunction
|
||||
{
|
||||
public:
|
||||
/**
|
||||
* Constructor.
|
||||
*
|
||||
* \param[in] pvtRegionIdx The pvt region index
|
||||
* \param[in] depth Depth nodes.
|
||||
* \param[in] pbub Dew-point pressure at @c depth.
|
||||
*/
|
||||
PDVD(const int pvtRegionIdx,
|
||||
const std::vector<double>& depth,
|
||||
const std::vector<double>& pdew)
|
||||
: pvtRegionIdx_(pvtRegionIdx)
|
||||
, pdewVsDepth_(depth, pdew)
|
||||
{}
|
||||
|
||||
/**
|
||||
* Function call.
|
||||
*
|
||||
* \param[in] depth Depth at which to calculate RV
|
||||
* value.
|
||||
*
|
||||
* \param[in] cellPress Pressure in the cell
|
||||
*
|
||||
* \param[in] temp Temperature at which to calculate RV
|
||||
* value.
|
||||
*
|
||||
* \return Vaporized oil-gas ratio (RV) at depth @c
|
||||
* depth and pressure @c press.
|
||||
*/
|
||||
double operator()(const double depth,
|
||||
const double cellPress,
|
||||
const double temp,
|
||||
const double satOil = 0.0) const
|
||||
{
|
||||
double press = cellPress;
|
||||
if (satOil <= 0.0) {
|
||||
if (pdewVsDepth_.xMin() > depth)
|
||||
press = pdewVsDepth_.valueAt(0);
|
||||
else if (pdewVsDepth_.xMax() < depth)
|
||||
press = pdewVsDepth_.valueAt(pdewVsDepth_.numSamples() - 1);
|
||||
else
|
||||
press = pdewVsDepth_.eval(depth, /*extrapolate=*/false);
|
||||
}
|
||||
return satRv(std::min(press, cellPress), temp);
|
||||
}
|
||||
|
||||
private:
|
||||
typedef Opm::Tabulated1DFunction<double> PdewVsDepthFunc;
|
||||
|
||||
const int pvtRegionIdx_;
|
||||
PdewVsDepthFunc pdewVsDepth_;
|
||||
|
||||
double satRv(const double press, const double temp) const
|
||||
{
|
||||
return FluidSystem::gasPvt().saturatedOilVaporizationFactor(pvtRegionIdx_, temp, press);
|
||||
}
|
||||
};
|
||||
|
||||
|
||||
/**
|
||||
* Type that implements "vaporized oil-gas ratio"
|
||||
|
||||
@@ -45,6 +45,8 @@
|
||||
#include <opm/parser/eclipse/EclipseState/Tables/TableManager.hpp>
|
||||
#include <opm/parser/eclipse/EclipseState/Tables/RsvdTable.hpp>
|
||||
#include <opm/parser/eclipse/EclipseState/Tables/RvvdTable.hpp>
|
||||
#include <opm/parser/eclipse/EclipseState/Tables/PbvdTable.hpp>
|
||||
#include <opm/parser/eclipse/EclipseState/Tables/PdvdTable.hpp>
|
||||
#include <opm/common/OpmLog/OpmLog.hpp>
|
||||
#include <opm/common/data/SimulationDataContainer.hpp>
|
||||
|
||||
@@ -974,7 +976,6 @@ public:
|
||||
// Create Rs functions.
|
||||
rsFunc_.reserve(rec.size());
|
||||
if (FluidSystem::enableDissolvedGas()) {
|
||||
const Opm::TableContainer& rsvdTables = tables.getRsvdTables();
|
||||
for (size_t i = 0; i < rec.size(); ++i) {
|
||||
if (eqlmap.cells(i).empty()) {
|
||||
rsFunc_.push_back(std::shared_ptr<Miscibility::RsVD<FluidSystem>>());
|
||||
@@ -982,14 +983,26 @@ public:
|
||||
}
|
||||
const int pvtIdx = regionPvtIdx_[i];
|
||||
if (!rec[i].liveOilInitConstantRs()) {
|
||||
if (rsvdTables.size() <= 0) {
|
||||
OPM_THROW(std::runtime_error, "Cannot initialise: RSVD table not available.");
|
||||
const Opm::TableContainer& rsvdTables = tables.getRsvdTables();
|
||||
const Opm::TableContainer& pbvdTables = tables.getPbvdTables();
|
||||
if (rsvdTables.size() > 0) {
|
||||
|
||||
const Opm::RsvdTable& rsvdTable = rsvdTables.getTable<Opm::RsvdTable>(i);
|
||||
std::vector<double> depthColumn = rsvdTable.getColumn("DEPTH").vectorCopy();
|
||||
std::vector<double> rsColumn = rsvdTable.getColumn("RS").vectorCopy();
|
||||
rsFunc_.push_back(std::make_shared<Miscibility::RsVD<FluidSystem>>(pvtIdx,
|
||||
depthColumn, rsColumn));
|
||||
} else if (pbvdTables.size() > 0) {
|
||||
const Opm::PbvdTable& pbvdTable = pbvdTables.getTable<Opm::PbvdTable>(i);
|
||||
std::vector<double> depthColumn = pbvdTable.getColumn("DEPTH").vectorCopy();
|
||||
std::vector<double> pbubColumn = pbvdTable.getColumn("PBUB").vectorCopy();
|
||||
rsFunc_.push_back(std::make_shared<Miscibility::PBVD<FluidSystem>>(pvtIdx,
|
||||
depthColumn, pbubColumn));
|
||||
|
||||
} else {
|
||||
OPM_THROW(std::runtime_error, "Cannot initialise: RSVD or PBVD table not available.");
|
||||
}
|
||||
const Opm::RsvdTable& rsvdTable = rsvdTables.getTable<Opm::RsvdTable>(i);
|
||||
std::vector<double> depthColumn = rsvdTable.getColumn("DEPTH").vectorCopy();
|
||||
std::vector<double> rsColumn = rsvdTable.getColumn("RS").vectorCopy();
|
||||
rsFunc_.push_back(std::make_shared<Miscibility::RsVD<FluidSystem>>(pvtIdx,
|
||||
depthColumn, rsColumn));
|
||||
|
||||
}
|
||||
else {
|
||||
if (rec[i].gasOilContactDepth() != rec[i].datumDepth()) {
|
||||
@@ -1012,7 +1025,6 @@ public:
|
||||
|
||||
rvFunc_.reserve(rec.size());
|
||||
if (FluidSystem::enableVaporizedOil()) {
|
||||
const Opm::TableContainer& rvvdTables = tables.getRvvdTables();
|
||||
for (size_t i = 0; i < rec.size(); ++i) {
|
||||
if (eqlmap.cells(i).empty()) {
|
||||
rvFunc_.push_back(std::shared_ptr<Miscibility::RvVD<FluidSystem>>());
|
||||
@@ -1020,16 +1032,24 @@ public:
|
||||
}
|
||||
const int pvtIdx = regionPvtIdx_[i];
|
||||
if (!rec[i].wetGasInitConstantRv()) {
|
||||
if (rvvdTables.size() <= 0) {
|
||||
OPM_THROW(std::runtime_error, "Cannot initialise: RVVD table not available.");
|
||||
const Opm::TableContainer& rvvdTables = tables.getRvvdTables();
|
||||
const Opm::TableContainer& pdvdTables = tables.getPdvdTables();
|
||||
|
||||
if (rvvdTables.size() > 0) {
|
||||
const Opm::RvvdTable& rvvdTable = rvvdTables.getTable<Opm::RvvdTable>(i);
|
||||
std::vector<double> depthColumn = rvvdTable.getColumn("DEPTH").vectorCopy();
|
||||
std::vector<double> rvColumn = rvvdTable.getColumn("RV").vectorCopy();
|
||||
rvFunc_.push_back(std::make_shared<Miscibility::RvVD<FluidSystem>>(pvtIdx,
|
||||
depthColumn, rvColumn));
|
||||
} else if (pdvdTables.size() > 0) {
|
||||
const Opm::PdvdTable& pdvdTable = pdvdTables.getTable<Opm::PdvdTable>(i);
|
||||
std::vector<double> depthColumn = pdvdTable.getColumn("DEPTH").vectorCopy();
|
||||
std::vector<double> pdewColumn = pdvdTable.getColumn("PDEW").vectorCopy();
|
||||
rvFunc_.push_back(std::make_shared<Miscibility::PDVD<FluidSystem>>(pvtIdx,
|
||||
depthColumn, pdewColumn));
|
||||
} else {
|
||||
OPM_THROW(std::runtime_error, "Cannot initialise: RVVD or PDCD table not available.");
|
||||
}
|
||||
|
||||
const Opm::RvvdTable& rvvdTable = rvvdTables.getTable<Opm::RvvdTable>(i);
|
||||
std::vector<double> depthColumn = rvvdTable.getColumn("DEPTH").vectorCopy();
|
||||
std::vector<double> rvColumn = rvvdTable.getColumn("RV").vectorCopy();
|
||||
rvFunc_.push_back(std::make_shared<Miscibility::RvVD<FluidSystem>>(pvtIdx,
|
||||
depthColumn, rvColumn));
|
||||
|
||||
}
|
||||
else {
|
||||
if (rec[i].gasOilContactDepth() != rec[i].datumDepth()) {
|
||||
|
||||
Reference in New Issue
Block a user