2012-08-27 16:48:21 +02:00
|
|
|
/*
|
|
|
|
|
Copyright 2012 SINTEF ICT, Applied Mathematics.
|
|
|
|
|
|
|
|
|
|
This file is part of the Open Porous Media project (OPM).
|
|
|
|
|
|
|
|
|
|
OPM is free software: you can redistribute it and/or modify
|
|
|
|
|
it under the terms of the GNU General Public License as published by
|
|
|
|
|
the Free Software Foundation, either version 3 of the License, or
|
|
|
|
|
(at your option) any later version.
|
|
|
|
|
|
|
|
|
|
OPM is distributed in the hope that it will be useful,
|
|
|
|
|
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
|
|
|
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
|
|
|
|
GNU General Public License for more details.
|
|
|
|
|
|
|
|
|
|
You should have received a copy of the GNU General Public License
|
|
|
|
|
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
|
|
|
|
*/
|
|
|
|
|
|
|
|
|
|
|
2013-04-10 12:56:14 +02:00
|
|
|
#include "config.h"
|
2014-04-24 17:07:41 +02:00
|
|
|
#include <opm/core/props/pvt/PvtDead.hpp>
|
2016-01-20 14:12:04 +01:00
|
|
|
#include <opm/parser/eclipse/EclipseState/Tables/PvdoTable.hpp>
|
|
|
|
|
#include <opm/parser/eclipse/EclipseState/Tables/PvdgTable.hpp>
|
2012-08-27 16:48:21 +02:00
|
|
|
#include <algorithm>
|
|
|
|
|
|
|
|
|
|
// Extra includes for debug dumping of tables.
|
|
|
|
|
// #include <boost/lexical_cast.hpp>
|
|
|
|
|
// #include <string>
|
|
|
|
|
// #include <fstream>
|
|
|
|
|
|
|
|
|
|
namespace Opm
|
|
|
|
|
{
|
|
|
|
|
|
|
|
|
|
//------------------------------------------------------------------------
|
|
|
|
|
// Member functions
|
|
|
|
|
//-------------------------------------------------------------------------
|
2014-01-24 12:03:21 +01:00
|
|
|
/// Constructor
|
2015-09-21 18:22:09 +02:00
|
|
|
void PvtDead::initFromOil(const TableContainer& pvdoTables)
|
2014-01-24 12:03:21 +01:00
|
|
|
{
|
2014-09-17 12:44:40 +02:00
|
|
|
int numRegions = pvdoTables.size();
|
2014-04-24 17:07:41 +02:00
|
|
|
|
|
|
|
|
// resize the attributes of the object
|
|
|
|
|
b_.resize(numRegions);
|
|
|
|
|
viscosity_.resize(numRegions);
|
2014-09-22 14:21:33 +02:00
|
|
|
inverseBmu_.resize(numRegions);
|
2014-04-24 17:07:41 +02:00
|
|
|
|
|
|
|
|
for (int regionIdx = 0; regionIdx < numRegions; ++regionIdx) {
|
2015-09-21 18:22:09 +02:00
|
|
|
const Opm::PvdoTable& pvdoTable = pvdoTables.getTable<PvdoTable>(regionIdx);
|
2014-09-17 12:44:40 +02:00
|
|
|
|
2014-04-24 17:07:41 +02:00
|
|
|
// Copy data
|
2016-01-05 12:12:29 +01:00
|
|
|
const auto& press = pvdoTable.getColumn("P");
|
|
|
|
|
const auto& b_var = pvdoTable.getColumn("BO");
|
|
|
|
|
const auto& visc = pvdoTable.getColumn("MUO");
|
2014-04-24 17:07:41 +02:00
|
|
|
|
2015-09-02 13:55:46 +02:00
|
|
|
const int sz = b_var.size();
|
2014-09-22 14:21:33 +02:00
|
|
|
std::vector<double> inverseB(sz);
|
2014-04-24 17:07:41 +02:00
|
|
|
for (int i = 0; i < sz; ++i) {
|
2015-09-02 13:55:46 +02:00
|
|
|
inverseB[i] = 1.0 / b_var[i];
|
2014-04-24 17:07:41 +02:00
|
|
|
}
|
2014-09-12 13:39:16 +02:00
|
|
|
|
2014-09-22 14:21:33 +02:00
|
|
|
std::vector<double> inverseBmu(sz);
|
2014-09-12 13:39:16 +02:00
|
|
|
for (int i = 0; i < sz; ++i) {
|
2015-09-02 13:55:46 +02:00
|
|
|
inverseBmu[i] = 1.0 / (b_var[i] * visc[i]);
|
2014-09-12 13:39:16 +02:00
|
|
|
}
|
|
|
|
|
|
2014-09-22 14:21:33 +02:00
|
|
|
b_[regionIdx] = NonuniformTableLinear<double>(press, inverseB);
|
2014-04-24 17:07:41 +02:00
|
|
|
viscosity_[regionIdx] = NonuniformTableLinear<double>(press, visc);
|
2014-09-22 14:21:33 +02:00
|
|
|
inverseBmu_[regionIdx] = NonuniformTableLinear<double>(press, inverseBmu);
|
2014-01-24 12:03:21 +01:00
|
|
|
}
|
2014-03-03 09:17:11 +01:00
|
|
|
}
|
|
|
|
|
|
2014-04-24 17:07:41 +02:00
|
|
|
|
2015-09-21 18:22:09 +02:00
|
|
|
void PvtDead::initFromGas(const TableContainer& pvdgTables)
|
2014-03-03 09:17:11 +01:00
|
|
|
{
|
2014-09-17 12:44:40 +02:00
|
|
|
int numRegions = pvdgTables.size();
|
2014-04-24 17:07:41 +02:00
|
|
|
|
|
|
|
|
// resize the attributes of the object
|
|
|
|
|
b_.resize(numRegions);
|
|
|
|
|
viscosity_.resize(numRegions);
|
2014-09-22 14:21:33 +02:00
|
|
|
inverseBmu_.resize(numRegions);
|
2014-04-24 17:07:41 +02:00
|
|
|
|
|
|
|
|
for (int regionIdx = 0; regionIdx < numRegions; ++regionIdx) {
|
2015-09-21 18:22:09 +02:00
|
|
|
const Opm::PvdgTable& pvdgTable = pvdgTables.getTable<PvdgTable>(regionIdx);
|
2014-04-24 17:07:41 +02:00
|
|
|
|
|
|
|
|
// Copy data
|
2016-01-05 12:12:29 +01:00
|
|
|
const auto& press = pvdgTable.getColumn("P");
|
|
|
|
|
const auto& b = pvdgTable.getColumn("BG");
|
|
|
|
|
const auto& visc = pvdgTable.getColumn("MUG");
|
2014-04-24 17:07:41 +02:00
|
|
|
|
|
|
|
|
const int sz = b.size();
|
2014-09-22 14:21:33 +02:00
|
|
|
std::vector<double> inverseB(sz);
|
2014-04-24 17:07:41 +02:00
|
|
|
for (int i = 0; i < sz; ++i) {
|
2014-09-22 14:21:33 +02:00
|
|
|
inverseB[i] = 1.0 / b[i];
|
2014-04-24 17:07:41 +02:00
|
|
|
}
|
2014-09-12 13:39:16 +02:00
|
|
|
|
2014-09-22 14:21:33 +02:00
|
|
|
std::vector<double> inverseBmu(sz);
|
2014-09-12 13:39:16 +02:00
|
|
|
for (int i = 0; i < sz; ++i) {
|
2014-09-22 14:21:33 +02:00
|
|
|
inverseBmu[i] = 1.0 / (b[i] * visc[i]);
|
2014-09-12 13:39:16 +02:00
|
|
|
}
|
|
|
|
|
|
2014-09-22 14:21:33 +02:00
|
|
|
b_[regionIdx] = NonuniformTableLinear<double>(press, inverseB);
|
2014-04-24 17:07:41 +02:00
|
|
|
viscosity_[regionIdx] = NonuniformTableLinear<double>(press, visc);
|
2014-09-22 14:21:33 +02:00
|
|
|
inverseBmu_[regionIdx] = NonuniformTableLinear<double>(press, inverseBmu);
|
2014-03-03 09:17:11 +01:00
|
|
|
}
|
2012-08-27 16:48:21 +02:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
// Destructor
|
2014-04-24 17:07:41 +02:00
|
|
|
PvtDead::~PvtDead()
|
2012-08-27 16:48:21 +02:00
|
|
|
{
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
2014-04-24 17:07:41 +02:00
|
|
|
void PvtDead::mu(const int n,
|
2014-05-13 12:59:11 +02:00
|
|
|
const int* pvtTableIdx,
|
2014-11-20 12:15:01 +01:00
|
|
|
const double* p,
|
|
|
|
|
const double* /*T*/,
|
2012-08-27 16:48:21 +02:00
|
|
|
const double* /*z*/,
|
|
|
|
|
double* output_mu) const
|
|
|
|
|
{
|
|
|
|
|
// #pragma omp parallel for
|
|
|
|
|
for (int i = 0; i < n; ++i) {
|
2014-05-13 12:59:11 +02:00
|
|
|
int regionIdx = getTableIndex_(pvtTableIdx, i);
|
2014-09-14 17:33:43 +02:00
|
|
|
double tempInvB = b_[regionIdx](p[i]);
|
2014-09-22 14:21:33 +02:00
|
|
|
double tempInvBmu = inverseBmu_[regionIdx](p[i]);
|
|
|
|
|
output_mu[i] = tempInvB / tempInvBmu;
|
2012-08-27 16:48:21 +02:00
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
2014-04-24 17:07:41 +02:00
|
|
|
void PvtDead::mu(const int n,
|
2014-05-13 12:59:11 +02:00
|
|
|
const int* pvtTableIdx,
|
2014-11-20 12:15:01 +01:00
|
|
|
const double* p,
|
2015-03-17 12:40:05 +01:00
|
|
|
const double* /*T*/,
|
2013-05-22 12:57:01 +02:00
|
|
|
const double* /*r*/,
|
|
|
|
|
double* output_mu,
|
|
|
|
|
double* output_dmudp,
|
|
|
|
|
double* output_dmudr) const
|
|
|
|
|
{
|
|
|
|
|
// #pragma omp parallel for
|
|
|
|
|
for (int i = 0; i < n; ++i) {
|
2014-05-13 12:59:11 +02:00
|
|
|
int regionIdx = getTableIndex_(pvtTableIdx, i);
|
2014-09-14 17:33:43 +02:00
|
|
|
double tempInvB = b_[regionIdx](p[i]);
|
2014-09-22 14:21:33 +02:00
|
|
|
double tempInvBmu = inverseBmu_[regionIdx](p[i]);
|
|
|
|
|
output_mu[i] = tempInvB / tempInvBmu;
|
|
|
|
|
output_dmudp[i] = (tempInvBmu * b_[regionIdx].derivative(p[i])
|
|
|
|
|
- tempInvB * inverseBmu_[regionIdx].derivative(p[i])) / (tempInvBmu * tempInvBmu);
|
2013-05-22 12:57:01 +02:00
|
|
|
}
|
|
|
|
|
std::fill(output_dmudr, output_dmudr + n, 0.0);
|
|
|
|
|
|
|
|
|
|
}
|
|
|
|
|
|
2014-04-24 17:07:41 +02:00
|
|
|
void PvtDead::mu(const int n,
|
2014-05-13 12:59:11 +02:00
|
|
|
const int* pvtTableIdx,
|
2014-11-20 12:15:01 +01:00
|
|
|
const double* p,
|
2015-03-17 12:40:05 +01:00
|
|
|
const double* /*T*/,
|
2013-11-28 10:29:24 +01:00
|
|
|
const double* /*r*/,
|
2013-12-03 16:12:52 +01:00
|
|
|
const PhasePresence* /*cond*/,
|
2013-11-28 10:29:24 +01:00
|
|
|
double* output_mu,
|
|
|
|
|
double* output_dmudp,
|
|
|
|
|
double* output_dmudr) const
|
|
|
|
|
{
|
|
|
|
|
// #pragma omp parallel for
|
|
|
|
|
for (int i = 0; i < n; ++i) {
|
2014-05-13 12:59:11 +02:00
|
|
|
int regionIdx = getTableIndex_(pvtTableIdx, i);
|
2014-09-12 13:39:16 +02:00
|
|
|
double tempInvB = b_[regionIdx](p[i]);
|
2014-09-22 14:21:33 +02:00
|
|
|
double tempInvBmu = inverseBmu_[regionIdx](p[i]);
|
|
|
|
|
output_mu[i] = tempInvB / tempInvBmu;
|
|
|
|
|
output_dmudp[i] = (tempInvBmu * b_[regionIdx].derivative(p[i])
|
2014-09-22 14:26:00 +02:00
|
|
|
- tempInvB * inverseBmu_[regionIdx].derivative(p[i]))
|
|
|
|
|
/ (tempInvBmu * tempInvBmu);
|
2013-11-28 10:29:24 +01:00
|
|
|
}
|
|
|
|
|
std::fill(output_dmudr, output_dmudr + n, 0.0);
|
|
|
|
|
|
|
|
|
|
}
|
|
|
|
|
|
2014-04-24 17:07:41 +02:00
|
|
|
void PvtDead::B(const int n,
|
2014-05-13 12:59:11 +02:00
|
|
|
const int* pvtTableIdx,
|
2014-11-20 12:15:01 +01:00
|
|
|
const double* p,
|
|
|
|
|
const double* /*T*/,
|
2012-08-27 16:48:21 +02:00
|
|
|
const double* /*z*/,
|
|
|
|
|
double* output_B) const
|
|
|
|
|
{
|
|
|
|
|
// #pragma omp parallel for
|
2013-05-22 12:57:01 +02:00
|
|
|
// B = 1/b
|
2012-08-27 16:48:21 +02:00
|
|
|
for (int i = 0; i < n; ++i) {
|
2014-05-13 12:59:11 +02:00
|
|
|
int regionIdx = getTableIndex_(pvtTableIdx, i);
|
2014-04-24 17:07:41 +02:00
|
|
|
output_B[i] = 1.0/b_[regionIdx](p[i]);
|
2012-08-27 16:48:21 +02:00
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
2014-04-24 17:07:41 +02:00
|
|
|
void PvtDead::dBdp(const int n,
|
2014-05-13 12:59:11 +02:00
|
|
|
const int* pvtTableIdx,
|
2014-11-20 12:15:01 +01:00
|
|
|
const double* p,
|
|
|
|
|
const double* T,
|
2012-08-27 16:48:21 +02:00
|
|
|
const double* /*z*/,
|
|
|
|
|
double* output_B,
|
|
|
|
|
double* output_dBdp) const
|
|
|
|
|
{
|
2014-11-20 12:15:01 +01:00
|
|
|
B(n, pvtTableIdx, p, T, 0, output_B);
|
2012-08-27 16:48:21 +02:00
|
|
|
// #pragma omp parallel for
|
|
|
|
|
for (int i = 0; i < n; ++i) {
|
2014-05-13 12:59:11 +02:00
|
|
|
int regionIdx = getTableIndex_(pvtTableIdx, i);
|
2012-08-27 16:48:21 +02:00
|
|
|
double Bg = output_B[i];
|
2014-04-24 17:07:41 +02:00
|
|
|
output_dBdp[i] = -Bg*Bg*b_[regionIdx].derivative(p[i]);
|
2012-08-27 16:48:21 +02:00
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
2014-04-24 17:07:41 +02:00
|
|
|
void PvtDead::b(const int n,
|
2014-05-13 12:59:11 +02:00
|
|
|
const int* pvtTableIdx,
|
2014-11-20 12:15:01 +01:00
|
|
|
const double* p,
|
|
|
|
|
const double* /*T*/,
|
2013-05-22 12:57:01 +02:00
|
|
|
const double* /*r*/,
|
|
|
|
|
double* output_b,
|
|
|
|
|
double* output_dbdp,
|
|
|
|
|
double* output_dbdr) const
|
|
|
|
|
|
|
|
|
|
{
|
|
|
|
|
// #pragma omp parallel for
|
|
|
|
|
for (int i = 0; i < n; ++i) {
|
2014-05-13 12:59:11 +02:00
|
|
|
int regionIdx = getTableIndex_(pvtTableIdx, i);
|
2014-04-24 17:07:41 +02:00
|
|
|
|
|
|
|
|
output_b[i] = b_[regionIdx](p[i]);
|
|
|
|
|
output_dbdp[i] = b_[regionIdx].derivative(p[i]);
|
2013-11-28 10:29:24 +01:00
|
|
|
|
|
|
|
|
}
|
|
|
|
|
std::fill(output_dbdr, output_dbdr + n, 0.0);
|
|
|
|
|
|
|
|
|
|
}
|
|
|
|
|
|
2014-04-24 17:07:41 +02:00
|
|
|
void PvtDead::b(const int n,
|
2014-05-13 12:59:11 +02:00
|
|
|
const int* pvtTableIdx,
|
2014-11-20 12:15:01 +01:00
|
|
|
const double* p,
|
|
|
|
|
const double* /*T*/,
|
2013-11-28 10:29:24 +01:00
|
|
|
const double* /*r*/,
|
2013-12-03 16:12:52 +01:00
|
|
|
const PhasePresence* /*cond*/,
|
2013-11-28 10:29:24 +01:00
|
|
|
double* output_b,
|
|
|
|
|
double* output_dbdp,
|
|
|
|
|
double* output_dbdr) const
|
|
|
|
|
|
|
|
|
|
{
|
|
|
|
|
// #pragma omp parallel for
|
|
|
|
|
for (int i = 0; i < n; ++i) {
|
2014-05-13 12:59:11 +02:00
|
|
|
int regionIdx = getTableIndex_(pvtTableIdx, i);
|
2014-04-24 17:07:41 +02:00
|
|
|
|
|
|
|
|
output_b[i] = b_[regionIdx](p[i]);
|
|
|
|
|
output_dbdp[i] = b_[regionIdx].derivative(p[i]);
|
2013-05-22 12:57:01 +02:00
|
|
|
|
|
|
|
|
}
|
|
|
|
|
std::fill(output_dbdr, output_dbdr + n, 0.0);
|
|
|
|
|
|
|
|
|
|
}
|
|
|
|
|
|
2014-04-24 17:07:41 +02:00
|
|
|
void PvtDead::rsSat(const int n,
|
2014-08-27 19:43:53 +02:00
|
|
|
const int* /*pvtTableIdx*/,
|
2013-05-22 12:57:01 +02:00
|
|
|
const double* /*p*/,
|
2013-12-13 15:43:36 +01:00
|
|
|
double* output_rsSat,
|
|
|
|
|
double* output_drsSatdp) const
|
2013-05-22 12:57:01 +02:00
|
|
|
{
|
2013-12-13 15:43:36 +01:00
|
|
|
std::fill(output_rsSat, output_rsSat + n, 0.0);
|
|
|
|
|
std::fill(output_drsSatdp, output_drsSatdp + n, 0.0);
|
|
|
|
|
}
|
|
|
|
|
|
2014-04-24 17:07:41 +02:00
|
|
|
void PvtDead::rvSat(const int n,
|
2014-08-27 19:43:53 +02:00
|
|
|
const int* /*pvtTableIdx*/,
|
2013-12-13 15:43:36 +01:00
|
|
|
const double* /*p*/,
|
|
|
|
|
double* output_rvSat,
|
|
|
|
|
double* output_drvSatdp) const
|
|
|
|
|
{
|
|
|
|
|
std::fill(output_rvSat, output_rvSat + n, 0.0);
|
|
|
|
|
std::fill(output_drvSatdp, output_drvSatdp + n, 0.0);
|
2013-05-22 12:57:01 +02:00
|
|
|
}
|
2012-08-27 16:48:21 +02:00
|
|
|
|
2014-04-24 17:07:41 +02:00
|
|
|
void PvtDead::R(const int n,
|
2014-08-27 19:43:53 +02:00
|
|
|
const int* /*pvtTableIdx*/,
|
2012-08-27 16:48:21 +02:00
|
|
|
const double* /*p*/,
|
|
|
|
|
const double* /*z*/,
|
|
|
|
|
double* output_R) const
|
|
|
|
|
{
|
|
|
|
|
std::fill(output_R, output_R + n, 0.0);
|
|
|
|
|
}
|
|
|
|
|
|
2014-04-24 17:07:41 +02:00
|
|
|
void PvtDead::dRdp(const int n,
|
2014-08-27 19:43:53 +02:00
|
|
|
const int* /*pvtTableIdx*/,
|
2012-08-27 16:48:21 +02:00
|
|
|
const double* /*p*/,
|
|
|
|
|
const double* /*z*/,
|
|
|
|
|
double* output_R,
|
|
|
|
|
double* output_dRdp) const
|
|
|
|
|
{
|
|
|
|
|
std::fill(output_R, output_R + n, 0.0);
|
|
|
|
|
std::fill(output_dRdp, output_dRdp + n, 0.0);
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
}
|