Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
/*
|
|
|
|
Copyright 2014 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/>.
|
|
|
|
*/
|
|
|
|
|
|
|
|
#ifndef OPM_INITSTATEEQUIL_HEADER_INCLUDED
|
|
|
|
#define OPM_INITSTATEEQUIL_HEADER_INCLUDED
|
|
|
|
|
2014-01-23 03:16:49 -06:00
|
|
|
#include <opm/core/io/eclipse/EclipseGridParser.hpp>
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
#include <opm/core/props/BlackoilPropertiesInterface.hpp>
|
|
|
|
#include <opm/core/props/BlackoilPhases.hpp>
|
|
|
|
#include <opm/core/utility/linearInterpolation.hpp>
|
2014-02-24 08:19:04 -06:00
|
|
|
#include <opm/core/utility/RegionMapping.hpp>
|
2014-02-19 06:42:07 -06:00
|
|
|
#include <opm/core/utility/RootFinders.hpp>
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
#include <opm/core/utility/Units.hpp>
|
|
|
|
|
|
|
|
#include <array>
|
|
|
|
#include <cassert>
|
2014-01-17 13:07:51 -06:00
|
|
|
#include <utility>
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
#include <vector>
|
|
|
|
|
2014-01-19 18:25:33 -06:00
|
|
|
/**
|
|
|
|
* \file
|
|
|
|
* Facilities for an ECLIPSE-style equilibration-based
|
|
|
|
* initialisation scheme (keyword 'EQUIL').
|
|
|
|
*/
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
struct UnstructuredGrid;
|
|
|
|
|
|
|
|
namespace Opm
|
|
|
|
{
|
2014-01-19 18:25:33 -06:00
|
|
|
/**
|
|
|
|
* Types and routines that collectively implement a basic
|
|
|
|
* ECLIPSE-style equilibration-based initialisation scheme.
|
|
|
|
*
|
|
|
|
* This namespace is intentionally nested to avoid name clashes
|
|
|
|
* with other parts of OPM.
|
|
|
|
*/
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
namespace equil {
|
|
|
|
template <class Props>
|
|
|
|
class DensityCalculator;
|
|
|
|
|
2014-01-19 18:25:33 -06:00
|
|
|
/**
|
|
|
|
* Facility for calculating phase densities based on the
|
|
|
|
* BlackoilPropertiesInterface.
|
|
|
|
*
|
|
|
|
* Implements the crucial <CODE>operator()(p,svol)</CODE>
|
|
|
|
* function that is expected by class EquilReg.
|
|
|
|
*/
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
template <>
|
|
|
|
class DensityCalculator< BlackoilPropertiesInterface > {
|
|
|
|
public:
|
2014-01-19 18:25:33 -06:00
|
|
|
/**
|
|
|
|
* Constructor.
|
|
|
|
*
|
|
|
|
* \param[in] props Implementation of the
|
|
|
|
* BlackoilPropertiesInterface.
|
|
|
|
*
|
|
|
|
* \param[in] c Single cell used as a representative cell
|
|
|
|
* in a PVT region.
|
|
|
|
*/
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
DensityCalculator(const BlackoilPropertiesInterface& props,
|
|
|
|
const int c)
|
|
|
|
: props_(props)
|
|
|
|
, c_(1, c)
|
|
|
|
{
|
|
|
|
}
|
|
|
|
|
2014-01-19 18:25:33 -06:00
|
|
|
/**
|
|
|
|
* Compute phase densities of all phases at phase point
|
|
|
|
* given by (pressure, surface volume) tuple.
|
|
|
|
*
|
|
|
|
* \param[in] p Fluid pressure.
|
|
|
|
*
|
|
|
|
* \param[in] z Surface volumes of all phases.
|
|
|
|
*
|
|
|
|
* \return Phase densities at phase point.
|
|
|
|
*/
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
std::vector<double>
|
|
|
|
operator()(const double p,
|
|
|
|
const std::vector<double>& z) const
|
|
|
|
{
|
|
|
|
const int np = props_.numPhases();
|
|
|
|
std::vector<double> A(np * np, 0);
|
|
|
|
|
|
|
|
assert (z.size() == std::vector<double>::size_type(np));
|
|
|
|
|
|
|
|
double* dAdp = 0;
|
|
|
|
props_.matrix(1, &p, &z[0], &c_[0], &A[0], dAdp);
|
|
|
|
|
|
|
|
std::vector<double> rho(np, 0.0);
|
|
|
|
props_.density(1, &A[0], &rho[0]);
|
|
|
|
|
|
|
|
return rho;
|
|
|
|
}
|
|
|
|
|
|
|
|
private:
|
|
|
|
const BlackoilPropertiesInterface& props_;
|
|
|
|
const std::vector<int> c_;
|
|
|
|
};
|
|
|
|
|
2014-01-19 18:25:33 -06:00
|
|
|
/**
|
|
|
|
* Types and routines relating to phase mixing in
|
|
|
|
* equilibration calculations.
|
|
|
|
*/
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
namespace miscibility {
|
2014-01-19 18:25:33 -06:00
|
|
|
/**
|
|
|
|
* Type that implements "no phase mixing" policy.
|
|
|
|
*/
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
struct NoMixing {
|
2014-01-19 18:25:33 -06:00
|
|
|
/**
|
|
|
|
* Function call.
|
|
|
|
*
|
|
|
|
* \param[in] depth Depth at which to calculate RS
|
|
|
|
* value.
|
|
|
|
*
|
|
|
|
* \param[in] press Pressure at which to calculate RS
|
|
|
|
* value.
|
|
|
|
*
|
|
|
|
* \return Dissolved gas-oil ratio (RS) at depth @c
|
|
|
|
* depth and pressure @c press. In "no mixing
|
|
|
|
* policy", this is identically zero.
|
|
|
|
*/
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
double
|
|
|
|
operator()(const double /* depth */,
|
|
|
|
const double /* press */) const
|
|
|
|
{
|
|
|
|
return 0.0;
|
|
|
|
}
|
|
|
|
};
|
|
|
|
|
2014-01-19 18:25:33 -06:00
|
|
|
/**
|
|
|
|
* Type that implements "dissolved gas-oil ratio"
|
|
|
|
* tabulated as a function of depth policy. Data
|
|
|
|
* typically taken from keyword 'RSVD'.
|
|
|
|
*/
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
class RsVD {
|
|
|
|
public:
|
2014-01-19 18:25:33 -06:00
|
|
|
/**
|
|
|
|
* Constructor.
|
|
|
|
*
|
|
|
|
* \param[in] depth Depth nodes.
|
|
|
|
* \param[in] rs Dissolved gas-oil ratio at @c depth.
|
|
|
|
*/
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
RsVD(const std::vector<double>& depth,
|
|
|
|
const std::vector<double>& rs)
|
|
|
|
: depth_(depth)
|
|
|
|
, rs_(rs)
|
|
|
|
{
|
|
|
|
}
|
|
|
|
|
2014-01-19 18:25:33 -06:00
|
|
|
/**
|
|
|
|
* Function call.
|
|
|
|
*
|
|
|
|
* \param[in] depth Depth at which to calculate RS
|
|
|
|
* value.
|
|
|
|
*
|
|
|
|
* \param[in] press Pressure at which to calculate RS
|
|
|
|
* value.
|
|
|
|
*
|
|
|
|
* \return Dissolved gas-oil ratio (RS) at depth @c
|
|
|
|
* depth and pressure @c press.
|
|
|
|
*/
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
double
|
|
|
|
operator()(const double depth,
|
|
|
|
const double /* press */) const
|
|
|
|
{
|
|
|
|
return linearInterpolation(depth_, rs_, depth);
|
|
|
|
}
|
|
|
|
|
|
|
|
private:
|
2014-01-19 18:25:33 -06:00
|
|
|
std::vector<double> depth_; /**< Depth nodes */
|
|
|
|
std::vector<double> rs_; /**< Dissolved gas-oil ratio */
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
};
|
2014-02-24 06:47:03 -06:00
|
|
|
|
|
|
|
|
|
|
|
/**
|
|
|
|
* Class that implements "dissolved gas-oil ratio" (Rs)
|
|
|
|
* as function of depth and pressure as follows:
|
|
|
|
*
|
|
|
|
* 1. The Rs at the gas-oil contact is equal to the
|
|
|
|
* saturated Rs value, Rs_sat_contact.
|
|
|
|
*
|
|
|
|
* 2. The Rs elsewhere is equal to Rs_sat_contact, but
|
|
|
|
* constrained to the saturated value as given by the
|
|
|
|
* local pressure.
|
|
|
|
*
|
|
|
|
* This should yield Rs-values that are constant below the
|
|
|
|
* contact, and decreasing above the contact.
|
|
|
|
*/
|
|
|
|
class RsSatAtContact {
|
|
|
|
public:
|
|
|
|
/**
|
|
|
|
* Constructor.
|
|
|
|
*
|
|
|
|
* \param[in] props property object
|
|
|
|
* \param[in] cell any cell in the pvt region
|
|
|
|
* \param[in] p_contact oil pressure at the contact
|
|
|
|
*/
|
|
|
|
RsSatAtContact(const BlackoilPropertiesInterface& props, const int cell, const double p_contact)
|
|
|
|
: props_(props), cell_(cell)
|
|
|
|
{
|
|
|
|
auto pu = props_.phaseUsage();
|
|
|
|
std::fill(z_, z_ + BlackoilPhases::MaxNumPhases, 0.0);
|
|
|
|
z_[pu.phase_pos[BlackoilPhases::Vapour]] = 1e100;
|
|
|
|
z_[pu.phase_pos[BlackoilPhases::Liquid]] = 1.0;
|
|
|
|
rs_sat_contact_ = satRs(p_contact);
|
|
|
|
}
|
|
|
|
|
|
|
|
/**
|
|
|
|
* Function call.
|
|
|
|
*
|
|
|
|
* \param[in] depth Depth at which to calculate RS
|
|
|
|
* value.
|
|
|
|
*
|
|
|
|
* \param[in] press Pressure 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 press) const
|
|
|
|
{
|
|
|
|
return std::max(satRs(press), rs_sat_contact_);
|
|
|
|
}
|
|
|
|
|
|
|
|
private:
|
|
|
|
const BlackoilPropertiesInterface& props_;
|
|
|
|
const int cell_;
|
|
|
|
double z_[BlackoilPhases::MaxNumPhases];
|
|
|
|
double rs_sat_contact_;
|
|
|
|
mutable double A_[BlackoilPhases::MaxNumPhases * BlackoilPhases::MaxNumPhases];
|
|
|
|
|
|
|
|
double satRs(const double press) const
|
|
|
|
{
|
|
|
|
props_.matrix(1, &press, z_, &cell_, A_, 0);
|
|
|
|
// Rs/Bo is in the gas row and oil column of A_.
|
|
|
|
// 1/Bo is in the oil row and column.
|
|
|
|
// Recall also that it is stored in column-major order.
|
|
|
|
const int opos = props_.phaseUsage().phase_pos[BlackoilPhases::Liquid];
|
|
|
|
const int gpos = props_.phaseUsage().phase_pos[BlackoilPhases::Vapour];
|
|
|
|
const int np = props_.numPhases();
|
|
|
|
return A_[np*opos + gpos] / A_[np*opos + opos];
|
|
|
|
}
|
|
|
|
};
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
} // namespace miscibility
|
|
|
|
|
2014-02-24 06:47:03 -06:00
|
|
|
|
2014-01-17 13:07:51 -06:00
|
|
|
|
2014-01-19 18:25:33 -06:00
|
|
|
/**
|
|
|
|
* Equilibration record.
|
|
|
|
*
|
|
|
|
* Layout and contents inspired by first six items of
|
|
|
|
* ECLIPSE's 'EQUIL' records. This is the minimum amount of
|
|
|
|
* input data needed to define phase pressures in an
|
|
|
|
* equilibration region.
|
|
|
|
*
|
|
|
|
* Data consists of three pairs of depth and pressure values:
|
|
|
|
* 1. main
|
|
|
|
* - @c depth Main datum depth.
|
|
|
|
* - @c press Pressure at datum depth.
|
|
|
|
*
|
|
|
|
* 2. woc
|
|
|
|
* - @c depth Depth of water-oil contact
|
|
|
|
* - @c press water-oil capillary pressure at water-oil contact.
|
|
|
|
* Capillary pressure defined as "P_oil - P_water".
|
|
|
|
*
|
|
|
|
* 3. goc
|
|
|
|
* - @c depth Depth of gas-oil contact
|
|
|
|
* - @c press Gas-oil capillary pressure at gas-oil contact.
|
|
|
|
* Capillary pressure defined as "P_gas - P_oil".
|
|
|
|
*/
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
struct EquilRecord {
|
|
|
|
struct {
|
|
|
|
double depth;
|
|
|
|
double press;
|
|
|
|
} main, woc, goc;
|
|
|
|
};
|
|
|
|
|
2014-01-19 18:25:33 -06:00
|
|
|
/**
|
|
|
|
* Aggregate information base of an equilibration region.
|
|
|
|
*
|
|
|
|
* Provides inquiry methods for retrieving depths of contacs
|
|
|
|
* and pressure values as well as a means of calculating fluid
|
|
|
|
* densities, dissolved gas-oil ratio and vapourised oil-gas
|
|
|
|
* ratios.
|
|
|
|
*
|
|
|
|
* \tparam DensCalc Type that provides access to a phase
|
|
|
|
* density calculation facility. Must implement an operator()
|
|
|
|
* declared as
|
|
|
|
* <CODE>
|
|
|
|
* std::vector<double>
|
|
|
|
* operator()(const double press,
|
|
|
|
* const std::vector<double>& svol )
|
|
|
|
* </CODE>
|
|
|
|
* that calculates the phase densities of all phases in @c
|
|
|
|
* svol at fluid pressure @c press.
|
|
|
|
*
|
|
|
|
* \tparam RS Type that provides access to a calculator for
|
|
|
|
* (initial) dissolved gas-oil ratios as a function of depth
|
|
|
|
* and (oil) pressure. Must implement an operator() declared
|
|
|
|
* as
|
|
|
|
* <CODE>
|
|
|
|
* double
|
|
|
|
* operator()(const double depth,
|
|
|
|
* const double press)
|
|
|
|
* </CODE>
|
|
|
|
* that calculates the dissolved gas-oil ratio at depth @c
|
|
|
|
* depth and (oil) pressure @c press.
|
|
|
|
*
|
|
|
|
* \tparam RV Type that provides access to a calculator for
|
|
|
|
* (initial) vapourised oil-gas ratios as a function of depth
|
|
|
|
* and (gas) pressure. Must implement an operator() declared
|
|
|
|
* as
|
|
|
|
* <CODE>
|
|
|
|
* double
|
|
|
|
* operator()(const double depth,
|
|
|
|
* const double press)
|
|
|
|
* </CODE>
|
|
|
|
* that calculates the vapourised oil-gas ratio at depth @c
|
|
|
|
* depth and (gas) pressure @c press.
|
|
|
|
*/
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
template <class DensCalc,
|
|
|
|
class RS = miscibility::NoMixing,
|
|
|
|
class RV = miscibility::NoMixing>
|
|
|
|
class EquilReg {
|
|
|
|
public:
|
2014-01-19 18:25:33 -06:00
|
|
|
/**
|
|
|
|
* Constructor.
|
|
|
|
*
|
|
|
|
* \param[in] rec Equilibration data of current region.
|
|
|
|
* \param[in] density Density calculator of current region.
|
|
|
|
* \param[in] rs Calculator of dissolved gas-oil ratio.
|
|
|
|
* \param[in] rv Calculator of vapourised oil-gas ratio.
|
|
|
|
* \param[in] pu Summary of current active phases.
|
|
|
|
*/
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
EquilReg(const EquilRecord& rec,
|
|
|
|
const DensCalc& density,
|
|
|
|
const RS& rs,
|
|
|
|
const RV& rv,
|
|
|
|
const PhaseUsage& pu)
|
|
|
|
: rec_ (rec)
|
|
|
|
, density_(density)
|
|
|
|
, rs_ (rs)
|
|
|
|
, rv_ (rv)
|
|
|
|
, pu_ (pu)
|
|
|
|
{
|
|
|
|
}
|
|
|
|
|
2014-01-19 18:25:33 -06:00
|
|
|
/**
|
|
|
|
* Type of density calculator.
|
|
|
|
*/
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
typedef DensCalc CalcDensity;
|
|
|
|
|
2014-01-19 18:25:33 -06:00
|
|
|
/**
|
|
|
|
* Type of dissolved gas-oil ratio calculator.
|
|
|
|
*/
|
|
|
|
typedef RS CalcDissolution;
|
|
|
|
|
|
|
|
/**
|
|
|
|
* Type of vapourised oil-gas ratio calculator.
|
|
|
|
*/
|
|
|
|
typedef RV CalcEvaporation;
|
|
|
|
|
|
|
|
/**
|
|
|
|
* Datum depth in current region
|
|
|
|
*/
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
double datum() const { return this->rec_.main.depth; }
|
2014-01-19 18:25:33 -06:00
|
|
|
|
|
|
|
/**
|
|
|
|
* Pressure at datum depth in current region.
|
|
|
|
*/
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
double pressure() const { return this->rec_.main.press; }
|
|
|
|
|
2014-01-19 18:25:33 -06:00
|
|
|
/**
|
|
|
|
* Depth of water-oil contact.
|
|
|
|
*/
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
double zwoc() const { return this->rec_.woc .depth; }
|
2014-01-19 18:25:33 -06:00
|
|
|
|
|
|
|
/**
|
|
|
|
* water-oil capillary pressure at water-oil contact.
|
|
|
|
*
|
|
|
|
* \return P_o - P_w at WOC.
|
|
|
|
*/
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
double pcow_woc() const { return this->rec_.woc .press; }
|
|
|
|
|
2014-01-19 18:25:33 -06:00
|
|
|
/**
|
|
|
|
* Depth of gas-oil contact.
|
|
|
|
*/
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
double zgoc() const { return this->rec_.goc .depth; }
|
2014-01-19 18:25:33 -06:00
|
|
|
|
|
|
|
/**
|
|
|
|
* Gas-oil capillary pressure at gas-oil contact.
|
|
|
|
*
|
|
|
|
* \return P_g - P_o at GOC.
|
|
|
|
*/
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
double pcgo_goc() const { return this->rec_.goc .press; }
|
|
|
|
|
2014-01-19 18:25:33 -06:00
|
|
|
/**
|
|
|
|
* Retrieve phase density calculator of current region.
|
|
|
|
*/
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
const CalcDensity&
|
|
|
|
densityCalculator() const { return this->density_; }
|
|
|
|
|
2014-01-19 18:25:33 -06:00
|
|
|
/**
|
|
|
|
* Retrieve dissolved gas-oil ratio calculator of current
|
|
|
|
* region.
|
|
|
|
*/
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
const CalcDissolution&
|
|
|
|
dissolutionCalculator() const { return this->rs_; }
|
|
|
|
|
2014-01-19 18:25:33 -06:00
|
|
|
/**
|
|
|
|
* Retrieve vapourised oil-gas ratio calculator of current
|
|
|
|
* region.
|
|
|
|
*/
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
const CalcEvaporation&
|
|
|
|
evaporationCalculator() const { return this->rv_; }
|
|
|
|
|
2014-01-19 18:25:33 -06:00
|
|
|
/**
|
|
|
|
* Retrieve active fluid phase summary.
|
|
|
|
*/
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
const PhaseUsage&
|
|
|
|
phaseUsage() const { return this->pu_; }
|
|
|
|
|
|
|
|
private:
|
2014-01-19 18:25:33 -06:00
|
|
|
EquilRecord rec_; /**< Equilibration data */
|
|
|
|
DensCalc density_; /**< Density calculator */
|
|
|
|
RS rs_; /**< RS calculator */
|
|
|
|
RV rv_; /**< RV calculator */
|
|
|
|
PhaseUsage pu_; /**< Active phase summary */
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
};
|
|
|
|
|
2014-02-19 06:42:07 -06:00
|
|
|
|
|
|
|
|
|
|
|
/// Functor for inverting capillary pressure function.
|
|
|
|
/// Function represented is
|
|
|
|
/// f(s) = pc(s) - target_pc
|
|
|
|
struct PcEq
|
|
|
|
{
|
|
|
|
PcEq(const BlackoilPropertiesInterface& props,
|
|
|
|
const int phase,
|
|
|
|
const int cell,
|
|
|
|
const double target_pc)
|
|
|
|
: props_(props),
|
|
|
|
phase_(phase),
|
|
|
|
cell_(cell),
|
|
|
|
target_pc_(target_pc)
|
|
|
|
{
|
|
|
|
std::fill(s_, s_ + BlackoilPhases::MaxNumPhases, 0.0);
|
|
|
|
std::fill(pc_, pc_ + BlackoilPhases::MaxNumPhases, 0.0);
|
|
|
|
}
|
|
|
|
double operator()(double s) const
|
|
|
|
{
|
|
|
|
s_[phase_] = s;
|
|
|
|
props_.capPress(1, s_, &cell_, pc_, 0);
|
|
|
|
return pc_[phase_] - target_pc_;
|
|
|
|
}
|
|
|
|
private:
|
|
|
|
const BlackoilPropertiesInterface& props_;
|
|
|
|
const int phase_;
|
|
|
|
const int cell_;
|
|
|
|
const int target_pc_;
|
|
|
|
mutable double s_[BlackoilPhases::MaxNumPhases];
|
|
|
|
mutable double pc_[BlackoilPhases::MaxNumPhases];
|
|
|
|
};
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
/// Compute saturation of some phase corresponding to a given
|
|
|
|
/// capillary pressure.
|
|
|
|
inline double satFromPc(const BlackoilPropertiesInterface& props,
|
|
|
|
const int phase,
|
|
|
|
const int cell,
|
|
|
|
const double target_pc,
|
|
|
|
const bool increasing = false)
|
|
|
|
{
|
|
|
|
// Find minimum and maximum saturations.
|
|
|
|
double sminarr[BlackoilPhases::MaxNumPhases];
|
|
|
|
double smaxarr[BlackoilPhases::MaxNumPhases];
|
|
|
|
props.satRange(1, &cell, sminarr, smaxarr);
|
|
|
|
const double s0 = increasing ? smaxarr[phase] : sminarr[phase];
|
|
|
|
const double s1 = increasing ? sminarr[phase] : smaxarr[phase];
|
|
|
|
|
|
|
|
// Create the equation f(s) = pc(s) - target_pc
|
|
|
|
const PcEq f(props, phase, cell, target_pc);
|
|
|
|
const double f0 = f(s0);
|
|
|
|
const double f1 = f(s1);
|
|
|
|
if (f0 <= 0.0) {
|
|
|
|
return s0;
|
|
|
|
} else if (f1 > 0.0) {
|
|
|
|
return s1;
|
|
|
|
} else {
|
|
|
|
const int max_iter = 30;
|
|
|
|
const double tol = 1e-6;
|
|
|
|
int iter_used = -1;
|
|
|
|
typedef RegulaFalsi<ThrowOnError> ScalarSolver;
|
|
|
|
const double sol = ScalarSolver::solve(f, std::min(s0, s1), std::max(s0, s1), max_iter, tol, iter_used);
|
|
|
|
return sol;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
/// Functor for inverting a sum of capillary pressure functions.
|
|
|
|
/// Function represented is
|
|
|
|
/// f(s) = pc1(s) + pc2(1 - s) - target_pc
|
|
|
|
struct PcEqSum
|
|
|
|
{
|
|
|
|
PcEqSum(const BlackoilPropertiesInterface& props,
|
|
|
|
const int phase1,
|
|
|
|
const int phase2,
|
|
|
|
const int cell,
|
|
|
|
const double target_pc)
|
|
|
|
: props_(props),
|
|
|
|
phase1_(phase1),
|
|
|
|
phase2_(phase2),
|
|
|
|
cell_(cell),
|
|
|
|
target_pc_(target_pc)
|
|
|
|
{
|
|
|
|
std::fill(s_, s_ + BlackoilPhases::MaxNumPhases, 0.0);
|
|
|
|
std::fill(pc_, pc_ + BlackoilPhases::MaxNumPhases, 0.0);
|
|
|
|
}
|
|
|
|
double operator()(double s) const
|
|
|
|
{
|
|
|
|
s_[phase1_] = s;
|
|
|
|
s_[phase2_] = 1.0 - s;
|
|
|
|
props_.capPress(1, s_, &cell_, pc_, 0);
|
|
|
|
return pc_[phase1_] + pc_[phase2_] - target_pc_;
|
|
|
|
}
|
|
|
|
private:
|
|
|
|
const BlackoilPropertiesInterface& props_;
|
|
|
|
const int phase1_;
|
|
|
|
const int phase2_;
|
|
|
|
const int cell_;
|
|
|
|
const int target_pc_;
|
|
|
|
mutable double s_[BlackoilPhases::MaxNumPhases];
|
|
|
|
mutable double pc_[BlackoilPhases::MaxNumPhases];
|
|
|
|
};
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
/// Compute saturation of some phase corresponding to a given
|
|
|
|
/// capillary pressure, where the capillary pressure function
|
|
|
|
/// is given as a sum of two other functions.
|
|
|
|
inline double satFromSumOfPcs(const BlackoilPropertiesInterface& props,
|
|
|
|
const int phase1,
|
|
|
|
const int phase2,
|
|
|
|
const int cell,
|
|
|
|
const double target_pc)
|
|
|
|
{
|
|
|
|
// Find minimum and maximum saturations.
|
|
|
|
double sminarr[BlackoilPhases::MaxNumPhases];
|
|
|
|
double smaxarr[BlackoilPhases::MaxNumPhases];
|
|
|
|
props.satRange(1, &cell, sminarr, smaxarr);
|
|
|
|
const double smin = sminarr[phase1];
|
|
|
|
const double smax = smaxarr[phase1];
|
|
|
|
|
|
|
|
// Create the equation f(s) = pc1(s) + pc2(1-s) - target_pc
|
|
|
|
const PcEqSum f(props, phase1, phase2, cell, target_pc);
|
2014-02-20 08:24:27 -06:00
|
|
|
const double f0 = f(smin);
|
|
|
|
const double f1 = f(smax);
|
|
|
|
if (f0 <= 0.0) {
|
|
|
|
return smin;
|
|
|
|
} else if (f1 > 0.0) {
|
|
|
|
return smax;
|
|
|
|
} else {
|
|
|
|
const int max_iter = 30;
|
|
|
|
const double tol = 1e-6;
|
|
|
|
int iter_used = -1;
|
|
|
|
typedef RegulaFalsi<ThrowOnError> ScalarSolver;
|
|
|
|
const double sol = ScalarSolver::solve(f, smin, smax, max_iter, tol, iter_used);
|
|
|
|
return sol;
|
|
|
|
}
|
2014-02-19 06:42:07 -06:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
2014-01-19 18:25:33 -06:00
|
|
|
/**
|
|
|
|
* Compute initial phase pressures by means of equilibration.
|
|
|
|
*
|
|
|
|
* This function uses the information contained in an
|
|
|
|
* equilibration record (i.e., depths and pressurs) as well as
|
|
|
|
* a density calculator and related data to vertically
|
|
|
|
* integrate the phase pressure ODE
|
|
|
|
* \f[
|
|
|
|
* \frac{\mathrm{d}p_{\alpha}}{\mathrm{d}z} =
|
|
|
|
* \rho_{\alpha}(z,p_{\alpha})\cdot g
|
|
|
|
* \f]
|
|
|
|
* in which \f$\rho_{\alpha}$ denotes the fluid density of
|
|
|
|
* fluid phase \f$\alpha\f$, \f$p_{\alpha}\f$ is the
|
|
|
|
* corresponding phase pressure, \f$z\f$ is the depth and
|
|
|
|
* \f$g\f$ is the acceleration due to gravity (assumed
|
|
|
|
* directed downwords, in the positive \f$z\f$ direction).
|
|
|
|
*
|
|
|
|
* \tparam Region Type of an equilibration region information
|
|
|
|
* base. Typically an instance of the EquilReg
|
|
|
|
* class template.
|
|
|
|
*
|
|
|
|
* \tparam CellRange Type of cell range that demarcates the
|
|
|
|
* cells pertaining to the current
|
2014-01-20 03:33:42 -06:00
|
|
|
* equilibration region. Must implement
|
|
|
|
* methods begin() and end() to bound the range
|
|
|
|
* as well as provide an inner type,
|
|
|
|
* const_iterator, to traverse the range.
|
2014-01-19 18:25:33 -06:00
|
|
|
*
|
|
|
|
* \param[in] G Grid.
|
|
|
|
* \param[in] reg Current equilibration region.
|
|
|
|
* \param[in] cells Range that spans the cells of the current
|
|
|
|
* equilibration region.
|
|
|
|
* \param[in] grav Acceleration of gravity.
|
|
|
|
*
|
|
|
|
* \return Phase pressures, one vector for each active phase,
|
|
|
|
* of pressure values in each cell in the current
|
|
|
|
* equilibration region.
|
|
|
|
*/
|
2014-01-17 10:43:27 -06:00
|
|
|
template <class Region, class CellRange>
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
std::vector< std::vector<double> >
|
|
|
|
phasePressures(const UnstructuredGrid& G,
|
|
|
|
const Region& reg,
|
2014-01-17 10:43:27 -06:00
|
|
|
const CellRange& cells,
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
const double grav = unit::gravity);
|
2014-01-23 03:16:49 -06:00
|
|
|
|
2014-02-19 06:42:07 -06:00
|
|
|
|
|
|
|
|
|
|
|
/**
|
|
|
|
* Compute initial phase saturations by means of equilibration.
|
|
|
|
*
|
|
|
|
* \tparam Region Type of an equilibration region information
|
|
|
|
* base. Typically an instance of the EquilReg
|
|
|
|
* class template.
|
|
|
|
*
|
|
|
|
* \tparam CellRange Type of cell range that demarcates the
|
|
|
|
* cells pertaining to the current
|
|
|
|
* equilibration region. Must implement
|
|
|
|
* methods begin() and end() to bound the range
|
|
|
|
* as well as provide an inner type,
|
|
|
|
* const_iterator, to traverse the range.
|
|
|
|
*
|
|
|
|
* \param[in] reg Current equilibration region.
|
|
|
|
* \param[in] cells Range that spans the cells of the current
|
|
|
|
* equilibration region.
|
|
|
|
* \param[in] props Property object, needed for capillary functions.
|
|
|
|
* \param[in] phase_pressures Phase pressures, one vector for each active phase,
|
|
|
|
* of pressure values in each cell in the current
|
|
|
|
* equilibration region.
|
|
|
|
* \return Phase saturations, one vector for each phase, each containing
|
|
|
|
* one saturation value per cell in the region.
|
|
|
|
*/
|
|
|
|
template <class Region, class CellRange>
|
|
|
|
std::vector< std::vector<double> >
|
|
|
|
phaseSaturations(const Region& reg,
|
|
|
|
const CellRange& cells,
|
|
|
|
const BlackoilPropertiesInterface& props,
|
|
|
|
const std::vector< std::vector<double> >& phase_pressures)
|
|
|
|
{
|
|
|
|
const double z0 = reg.datum();
|
|
|
|
const double zwoc = reg.zwoc ();
|
|
|
|
const double zgoc = reg.zgoc ();
|
|
|
|
if ((zgoc > z0) || (z0 > zwoc)) {
|
|
|
|
OPM_THROW(std::runtime_error, "Cannot initialise: the datum depth must be in the oil zone.");
|
|
|
|
}
|
|
|
|
if (!reg.phaseUsage().phase_used[BlackoilPhases::Liquid]) {
|
|
|
|
OPM_THROW(std::runtime_error, "Cannot initialise: not handling water-gas cases.");
|
|
|
|
}
|
|
|
|
|
|
|
|
std::vector< std::vector<double> > phase_saturations = phase_pressures; // Just to get the right size.
|
2014-02-20 08:24:27 -06:00
|
|
|
double smin[BlackoilPhases::MaxNumPhases] = { 0.0 };
|
|
|
|
double smax[BlackoilPhases::MaxNumPhases] = { 0.0 };
|
2014-02-19 06:42:07 -06:00
|
|
|
|
2014-02-20 08:24:27 -06:00
|
|
|
const bool water = reg.phaseUsage().phase_used[BlackoilPhases::Aqua];
|
|
|
|
const bool gas = reg.phaseUsage().phase_used[BlackoilPhases::Vapour];
|
2014-02-19 06:42:07 -06:00
|
|
|
const int oilpos = reg.phaseUsage().phase_pos[BlackoilPhases::Liquid];
|
2014-02-20 08:24:27 -06:00
|
|
|
const int waterpos = reg.phaseUsage().phase_pos[BlackoilPhases::Aqua];
|
|
|
|
const int gaspos = reg.phaseUsage().phase_pos[BlackoilPhases::Vapour];
|
2014-02-19 06:42:07 -06:00
|
|
|
std::vector<double>::size_type local_index = 0;
|
|
|
|
for (typename CellRange::const_iterator ci = cells.begin(); ci != cells.end(); ++ci, ++local_index) {
|
|
|
|
const int cell = *ci;
|
2014-02-20 08:24:27 -06:00
|
|
|
props.satRange(1, &cell, smin, smax);
|
2014-02-19 06:42:07 -06:00
|
|
|
// Find saturations from pressure differences by
|
|
|
|
// inverting capillary pressure functions.
|
|
|
|
double sw = 0.0;
|
2014-02-20 08:24:27 -06:00
|
|
|
if (water) {
|
2014-02-19 06:42:07 -06:00
|
|
|
const double pcov = phase_pressures[oilpos][local_index] - phase_pressures[waterpos][local_index];
|
|
|
|
sw = satFromPc(props, waterpos, cell, pcov);
|
|
|
|
phase_saturations[waterpos][local_index] = sw;
|
|
|
|
}
|
|
|
|
double sg = 0.0;
|
2014-02-20 08:24:27 -06:00
|
|
|
if (gas) {
|
2014-02-19 06:42:07 -06:00
|
|
|
// Note that pcog is defined to be (pg - po), not (po - pg).
|
|
|
|
const double pcog = phase_pressures[gaspos][local_index] - phase_pressures[oilpos][local_index];
|
|
|
|
const double increasing = true; // pcog(sg) expected to be increasing function
|
|
|
|
sg = satFromPc(props, gaspos, cell, pcog, increasing);
|
|
|
|
phase_saturations[gaspos][local_index] = sg;
|
|
|
|
}
|
2014-02-21 07:47:14 -06:00
|
|
|
if (gas && water && (sg + sw > 1.0)) {
|
2014-02-19 06:42:07 -06:00
|
|
|
// Overlapping gas-oil and oil-water transition
|
2014-02-21 07:47:14 -06:00
|
|
|
// zones can lead to unphysical saturations when
|
|
|
|
// treated as above. Must recalculate using gas-water
|
2014-02-19 06:42:07 -06:00
|
|
|
// capillary pressure.
|
|
|
|
const double pcgw = phase_pressures[gaspos][local_index] - phase_pressures[waterpos][local_index];
|
|
|
|
sw = satFromSumOfPcs(props, waterpos, gaspos, cell, pcgw);
|
|
|
|
sg = 1.0 - sw;
|
|
|
|
phase_saturations[waterpos][local_index] = sw;
|
|
|
|
phase_saturations[gaspos][local_index] = sg;
|
|
|
|
}
|
|
|
|
phase_saturations[oilpos][local_index] = 1.0 - sw - sg;
|
|
|
|
}
|
2014-02-20 08:24:27 -06:00
|
|
|
return phase_saturations;
|
2014-02-19 06:42:07 -06:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
2014-01-23 03:16:49 -06:00
|
|
|
namespace DeckDependent {
|
|
|
|
inline
|
|
|
|
std::vector<EquilRecord>
|
|
|
|
getEquil(const EclipseGridParser& deck)
|
|
|
|
{
|
|
|
|
if (deck.hasField("EQUIL")) {
|
|
|
|
const EQUIL& eql = deck.getEQUIL();
|
|
|
|
|
|
|
|
typedef std::vector<EquilLine>::size_type sz_t;
|
|
|
|
const sz_t nrec = eql.equil.size();
|
|
|
|
|
|
|
|
std::vector<EquilRecord> ret;
|
|
|
|
ret.reserve(nrec);
|
|
|
|
for (sz_t r = 0; r < nrec; ++r) {
|
|
|
|
const EquilLine& rec = eql.equil[r];
|
|
|
|
|
|
|
|
EquilRecord record =
|
|
|
|
{
|
|
|
|
{ rec.datum_depth_ ,
|
|
|
|
rec.datum_depth_pressure_ }
|
|
|
|
,
|
|
|
|
{ rec.water_oil_contact_depth_ ,
|
|
|
|
rec.oil_water_cap_pressure_ }
|
|
|
|
,
|
|
|
|
{ rec.gas_oil_contact_depth_ ,
|
|
|
|
rec.gas_oil_cap_pressure_ }
|
|
|
|
};
|
|
|
|
|
|
|
|
ret.push_back(record);
|
|
|
|
}
|
|
|
|
|
|
|
|
return ret;
|
|
|
|
}
|
|
|
|
else {
|
|
|
|
OPM_THROW(std::domain_error,
|
|
|
|
"Deck does not provide equilibration data.");
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
inline
|
|
|
|
std::vector<int>
|
|
|
|
equilnum(const EclipseGridParser& deck,
|
|
|
|
const UnstructuredGrid& G )
|
|
|
|
{
|
|
|
|
std::vector<int> eqlnum;
|
|
|
|
if (deck.hasField("EQLNUM")) {
|
|
|
|
eqlnum = deck.getIntegerValue("EQLNUM");
|
|
|
|
}
|
|
|
|
else {
|
|
|
|
// No explicit equilibration region.
|
|
|
|
// All cells in region zero.
|
|
|
|
eqlnum.assign(G.number_of_cells, 0);
|
|
|
|
}
|
|
|
|
|
|
|
|
return eqlnum;
|
|
|
|
}
|
|
|
|
|
|
|
|
template <class InputDeck>
|
2014-02-21 01:52:25 -06:00
|
|
|
class PhasePressureSaturationComputer;
|
2014-01-23 03:16:49 -06:00
|
|
|
|
|
|
|
template <>
|
2014-02-21 01:52:25 -06:00
|
|
|
class PhasePressureSaturationComputer<Opm::EclipseGridParser> {
|
2014-01-23 03:16:49 -06:00
|
|
|
public:
|
2014-02-21 01:52:25 -06:00
|
|
|
PhasePressureSaturationComputer(const BlackoilPropertiesInterface& props,
|
|
|
|
const EclipseGridParser& deck ,
|
|
|
|
const UnstructuredGrid& G ,
|
|
|
|
const double grav = unit::gravity)
|
2014-01-23 03:16:49 -06:00
|
|
|
: pp_(props.numPhases(),
|
2014-02-21 01:52:25 -06:00
|
|
|
std::vector<double>(G.number_of_cells)),
|
|
|
|
sat_(props.numPhases(),
|
2014-01-23 03:16:49 -06:00
|
|
|
std::vector<double>(G.number_of_cells))
|
|
|
|
{
|
|
|
|
const std::vector<EquilRecord> rec = getEquil(deck);
|
|
|
|
const RegionMapping<> eqlmap(equilnum(deck, G));
|
|
|
|
|
2014-02-21 01:52:25 -06:00
|
|
|
calcPressII(eqlmap, rec, props, G, grav);
|
|
|
|
calcSat(eqlmap, rec, props, G, grav);
|
2014-01-23 03:16:49 -06:00
|
|
|
}
|
|
|
|
|
|
|
|
typedef std::vector<double> PVal;
|
|
|
|
typedef std::vector<PVal> PPress;
|
|
|
|
|
|
|
|
const PPress& press() const { return pp_; }
|
2014-02-21 01:52:25 -06:00
|
|
|
const PPress& saturation() const { return sat_; }
|
2014-01-23 03:16:49 -06:00
|
|
|
|
|
|
|
private:
|
|
|
|
typedef DensityCalculator<BlackoilPropertiesInterface> RhoCalc;
|
|
|
|
typedef EquilReg<RhoCalc> EqReg;
|
|
|
|
|
|
|
|
PPress pp_;
|
2014-02-21 01:52:25 -06:00
|
|
|
PPress sat_;
|
2014-01-23 03:16:49 -06:00
|
|
|
|
|
|
|
template <class RMap>
|
|
|
|
void
|
2014-02-21 01:52:25 -06:00
|
|
|
calcPressII(const RMap& reg ,
|
|
|
|
const std::vector< EquilRecord >& rec ,
|
|
|
|
const Opm::BlackoilPropertiesInterface& props,
|
|
|
|
const UnstructuredGrid& G ,
|
|
|
|
const double grav)
|
2014-01-23 03:16:49 -06:00
|
|
|
{
|
|
|
|
typedef miscibility::NoMixing NoMix;
|
|
|
|
|
|
|
|
for (typename RMap::RegionId
|
|
|
|
r = 0, nr = reg.numRegions();
|
|
|
|
r < nr; ++r)
|
|
|
|
{
|
|
|
|
const typename RMap::CellRange cells = reg.cells(r);
|
|
|
|
|
|
|
|
const int repcell = *cells.begin();
|
|
|
|
const RhoCalc calc(props, repcell);
|
|
|
|
|
|
|
|
const EqReg eqreg(rec[r], calc, NoMix(), NoMix(),
|
|
|
|
props.phaseUsage());
|
|
|
|
|
2014-02-19 06:38:21 -06:00
|
|
|
const PPress& res = phasePressures(G, eqreg, cells, grav);
|
2014-01-23 03:16:49 -06:00
|
|
|
|
|
|
|
for (int p = 0, np = props.numPhases(); p < np; ++p) {
|
|
|
|
PVal& d = pp_[p];
|
|
|
|
PVal::const_iterator s = res[p].begin();
|
|
|
|
for (typename RMap::CellRange::const_iterator
|
|
|
|
c = cells.begin(),
|
|
|
|
e = cells.end();
|
|
|
|
c != e; ++c, ++s)
|
|
|
|
{
|
|
|
|
d[*c] = *s;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
2014-02-21 01:52:25 -06:00
|
|
|
|
|
|
|
template <class RMap>
|
|
|
|
void
|
|
|
|
calcSat(const RMap& reg ,
|
|
|
|
const std::vector< EquilRecord >& rec ,
|
|
|
|
const Opm::BlackoilPropertiesInterface& props,
|
|
|
|
const UnstructuredGrid& G ,
|
|
|
|
const double grav)
|
|
|
|
{
|
|
|
|
typedef miscibility::NoMixing NoMix;
|
|
|
|
|
|
|
|
for (typename RMap::RegionId
|
|
|
|
r = 0, nr = reg.numRegions();
|
|
|
|
r < nr; ++r)
|
|
|
|
{
|
|
|
|
const typename RMap::CellRange cells = reg.cells(r);
|
|
|
|
|
|
|
|
const int repcell = *cells.begin();
|
|
|
|
const RhoCalc calc(props, repcell);
|
|
|
|
|
|
|
|
const EqReg eqreg(rec[r], calc, NoMix(), NoMix(),
|
|
|
|
props.phaseUsage());
|
|
|
|
|
|
|
|
const PPress press = phasePressures(G, eqreg, cells, grav);
|
|
|
|
const PPress sat = phaseSaturations(eqreg, cells, props, press);
|
|
|
|
|
|
|
|
for (int p = 0, np = props.numPhases(); p < np; ++p) {
|
|
|
|
PVal& d = pp_[p];
|
|
|
|
PVal::const_iterator s = press[p].begin();
|
|
|
|
for (typename RMap::CellRange::const_iterator
|
|
|
|
c = cells.begin(),
|
|
|
|
e = cells.end();
|
|
|
|
c != e; ++c, ++s)
|
|
|
|
{
|
|
|
|
d[*c] = *s;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
for (int p = 0, np = props.numPhases(); p < np; ++p) {
|
|
|
|
PVal& d = sat_[p];
|
|
|
|
PVal::const_iterator s = sat[p].begin();
|
|
|
|
for (typename RMap::CellRange::const_iterator
|
|
|
|
c = cells.begin(),
|
|
|
|
e = cells.end();
|
|
|
|
c != e; ++c, ++s)
|
|
|
|
{
|
|
|
|
d[*c] = *s;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
|
2014-01-23 03:16:49 -06:00
|
|
|
};
|
|
|
|
} // namespace DeckDependent
|
Add basic equilibration facility
This commit adds a simple facility for calculating initial phase
pressures assuming stationary conditions, a known reference pressure
in the oil zone as well as the depth and capillary pressures at the
water-oil and gas-oil contacts.
Function 'Opm::equil::phasePressures()' uses a simple ODE/IVP-based
approach, solved using the traditional RK4 method with constant step
sizes, to derive the required pressure values. Specifically, we
solve the ODE
dp/dz = rho(z,p) * g
with 'z' represening depth, 'p' being a phase pressure and 'rho' the
associate phase density. Finally, 'g' is the acceleration of
gravity. We assume that we can calculate phase densities, e.g.,
from table look-up. This assumption holds in the case of an ECLIPSE
input deck.
Using RK4 with constant step sizes is a limitation of this
implementation. This, basically, assumes that the phase densities
varies only smoothly with depth and pressure (at reservoir
conditions).
2014-01-14 13:37:28 -06:00
|
|
|
} // namespace equil
|
|
|
|
} // namespace Opm
|
|
|
|
|
|
|
|
#include <opm/core/simulator/initStateEquil_impl.hpp>
|
|
|
|
|
|
|
|
#endif // OPM_INITSTATEEQUIL_HEADER_INCLUDED
|