Add Support Infrastructure for WELPI Feature
This commit adds logic implementing the static parts of the WELPI keyword. We internalize the keyword data, record appropriate events and provide hooks for dynamically adjusting the per-connection transmissibility factor (Connection::CF()) when those events occur. We implement support at three levels - WellConnections: Add new public member functions prepareWellPIScaling and applyWellPIScaling which, respectively, creates bookkeeping data to track those connections which are subject to CF scaling and actually applies that CF scaling. - Well: Add new data member 'productivity_index' which holds the 'WELPI' data value from the input keyword (converted to SI) and new member functions updateWellProductivityIndex and applyWellProdIndexScaling. The first follows the 'update*' protocol (return 'true' if state change) and assigns new values to 'productivity_index' while the second uses the stored PI value and a dynamically calculated effective PI value to rescale the pertinent connections' CF value. - Schedule: Add new member function handleWELPI which internalizes the WELPI keyword and its data and records WELPI events for subsequent playback in the simulator layer. Also add a set of unit tests to exercise the new features at all levels.
This commit is contained in:
parent
de5e3d90cd
commit
fa7d8bc28c
@ -101,6 +101,11 @@ namespace Opm
|
||||
*/
|
||||
GROUP_PRODUCTION_UPDATE = (1 << 14),
|
||||
GROUP_INJECTION_UPDATE = (1 << 15),
|
||||
|
||||
/*
|
||||
* New explicit well productivity/injectivity assignment.
|
||||
*/
|
||||
WELL_PRODUCTIVITY_INDEX = (1 << 16),
|
||||
};
|
||||
}
|
||||
|
||||
|
@ -521,6 +521,7 @@ namespace Opm
|
||||
void handleWECON (const HandlerContext&, const ParseContext&, ErrorGuard&);
|
||||
void handleWEFAC (const HandlerContext&, const ParseContext&, ErrorGuard&);
|
||||
void handleWELOPEN (const HandlerContext&, const ParseContext&, ErrorGuard&);
|
||||
void handleWELPI (const HandlerContext&, const ParseContext&, ErrorGuard&);
|
||||
void handleWELSEGS (const HandlerContext&, const ParseContext&, ErrorGuard&);
|
||||
void handleWELSPECS (const HandlerContext&, const ParseContext&, ErrorGuard&);
|
||||
void handleWELTARG (const HandlerContext&, const ParseContext&, ErrorGuard&);
|
||||
|
@ -23,6 +23,7 @@
|
||||
|
||||
#include <cstddef>
|
||||
#include <iosfwd>
|
||||
#include <limits>
|
||||
#include <map>
|
||||
#include <memory>
|
||||
#include <string>
|
||||
@ -544,6 +545,7 @@ public:
|
||||
bool updateEconLimits(std::shared_ptr<WellEconProductionLimits> econ_limits);
|
||||
bool updateProduction(std::shared_ptr<WellProductionProperties> production);
|
||||
bool updateInjection(std::shared_ptr<WellInjectionProperties> injection);
|
||||
bool updateWellProductivityIndex(const double prodIndex);
|
||||
bool updateWSEGSICD(const std::vector<std::pair<int, SICD> >& sicd_pairs);
|
||||
bool updateWSEGVALV(const std::vector<std::pair<int, Valve> >& valve_pairs);
|
||||
|
||||
@ -568,6 +570,7 @@ public:
|
||||
bool cmp_structure(const Well& other) const;
|
||||
bool operator==(const Well& data) const;
|
||||
void setInsertIndex(std::size_t index);
|
||||
void applyWellProdIndexScaling(const double currentEffectivePI);
|
||||
|
||||
template<class Serializer>
|
||||
void serializeOp(Serializer& serializer)
|
||||
@ -593,6 +596,7 @@ public:
|
||||
serializer(solvent_fraction);
|
||||
serializer(has_produced);
|
||||
serializer(prediction_mode);
|
||||
serializer(productivity_index);
|
||||
serializer(econ_limits);
|
||||
serializer(foam_properties);
|
||||
serializer(polymer_properties);
|
||||
@ -629,14 +633,14 @@ private:
|
||||
double solvent_fraction;
|
||||
bool has_produced = false;
|
||||
bool prediction_mode = true;
|
||||
|
||||
double productivity_index{ std::numeric_limits<double>::lowest() };
|
||||
|
||||
std::shared_ptr<WellEconProductionLimits> econ_limits;
|
||||
std::shared_ptr<WellFoamProperties> foam_properties;
|
||||
std::shared_ptr<WellPolymerProperties> polymer_properties;
|
||||
std::shared_ptr<WellBrineProperties> brine_properties;
|
||||
std::shared_ptr<WellTracerProperties> tracer_properties;
|
||||
std::shared_ptr<WellConnections> connections; // The WellConnections object can not be const because of the filterConnections method - would be beneficial to rewrite to enable const
|
||||
std::shared_ptr<WellConnections> connections; // The WellConnections object cannot be const because of WELPI and the filterConnections method
|
||||
std::shared_ptr<WellProductionProperties> production;
|
||||
std::shared_ptr<WellInjectionProperties> injection;
|
||||
std::shared_ptr<WellSegments> segments;
|
||||
|
@ -25,9 +25,13 @@
|
||||
|
||||
#include <opm/common/utility/ActiveGridCells.hpp>
|
||||
|
||||
#include <cstddef>
|
||||
#include <vector>
|
||||
|
||||
#include <stddef.h>
|
||||
|
||||
namespace Opm {
|
||||
class DeckRecord;
|
||||
class EclipseGrid;
|
||||
class FieldPropsManager;
|
||||
class WellConnections {
|
||||
@ -95,6 +99,21 @@ namespace Opm {
|
||||
Connection::Order ordering() const { return this->m_ordering; }
|
||||
std::vector<const Connection *> output(const EclipseGrid& grid) const;
|
||||
|
||||
/// Activate or reactivate WELPI scaling for this connection set.
|
||||
///
|
||||
/// Following this call, any WELPI-based scaling will apply to all
|
||||
/// connections whose properties are not reset in COMPDAT.
|
||||
///
|
||||
/// Returns whether or not this call to prepareWellPIScaling() is
|
||||
/// a state change (e.g., no WELPI to active WELPI or WELPI for
|
||||
/// some connections to WELPI for all connections).
|
||||
bool prepareWellPIScaling();
|
||||
|
||||
/// Scale pertinent connections' CF value by supplied value. Scaling
|
||||
/// factor typically derived from 'WELPI' input keyword and a dynamic
|
||||
/// productivity index calculation.
|
||||
void applyWellPIScaling(const double scaleFactor);
|
||||
|
||||
template<class Serializer>
|
||||
void serializeOp(Serializer& serializer)
|
||||
{
|
||||
@ -102,6 +121,8 @@ namespace Opm {
|
||||
serializer(headI);
|
||||
serializer(headJ);
|
||||
serializer.vector(m_connections);
|
||||
serializer(m_hasWellPIAdjustment);
|
||||
serializer(m_wellPIConnections);
|
||||
}
|
||||
|
||||
private:
|
||||
@ -133,9 +154,23 @@ namespace Opm {
|
||||
void orderTRACK();
|
||||
void orderMSW();
|
||||
|
||||
// Exclude specific connection from WELPI CF scaling. No action unless
|
||||
// this connection set has been prepared for WELPI.
|
||||
void excludeFromWellPI(const std::size_t connID);
|
||||
|
||||
Connection::Order m_ordering = Connection::Order::TRACK;
|
||||
int headI, headJ;
|
||||
std::vector< Connection > m_connections;
|
||||
|
||||
// Backing data for 'WELPI'.
|
||||
// 1. No adjustment if this set of connections has not been prepared
|
||||
// for WELPI (m_hasWellPIAdjustment == false, default value).
|
||||
//
|
||||
// 2. Otherwise, scale Connection::CF() by supplied scaling factor
|
||||
// for those connections that are marked in m_wellPIConnections.
|
||||
// Apply scaling to all connections if m_wellPIConnections.empty().
|
||||
bool m_hasWellPIAdjustment{false};
|
||||
std::vector<bool> m_wellPIConnections{};
|
||||
};
|
||||
}
|
||||
|
||||
|
@ -51,6 +51,7 @@
|
||||
#include <opm/parser/eclipse/Parser/ParserKeywords/W.hpp>
|
||||
|
||||
#include <opm/parser/eclipse/EclipseState/EclipseState.hpp>
|
||||
#include <opm/parser/eclipse/EclipseState/Runspec.hpp>
|
||||
#include <opm/parser/eclipse/EclipseState/Schedule/Action/ActionX.hpp>
|
||||
#include <opm/parser/eclipse/EclipseState/Schedule/Action/ActionResult.hpp>
|
||||
#include <opm/parser/eclipse/EclipseState/Schedule/DynamicState.hpp>
|
||||
@ -1065,6 +1066,48 @@ namespace {
|
||||
return applyWELOPEN(handlerContext.keyword, handlerContext.currentStep, parseContext, errors);
|
||||
}
|
||||
|
||||
void Schedule::handleWELPI(const HandlerContext& handlerContext, const ParseContext& parseContext, ErrorGuard& errors) {
|
||||
// Keyword structure
|
||||
//
|
||||
// WELPI
|
||||
// W1 123.45 /
|
||||
// W2* 456.78 /
|
||||
// *P 111.222 /
|
||||
// **X* 333.444 /
|
||||
// /
|
||||
//
|
||||
// Interpretation of productivity index (item 2) depends on well's preferred phase.
|
||||
using WELL_NAME = ParserKeywords::WELPI::WELL_NAME;
|
||||
using PI = ParserKeywords::WELPI::STEADY_STATE_PRODUCTIVITY_OR_INJECTIVITY_INDEX_VALUE;
|
||||
|
||||
const auto& usys = handlerContext.section.unitSystem();
|
||||
const auto gasPI = UnitSystem::measure::gas_productivity_index;
|
||||
const auto liqPI = UnitSystem::measure::liquid_productivity_index;
|
||||
|
||||
for (const auto& record : handlerContext.keyword) {
|
||||
const auto well_names = this->wellNames(record.getItem<WELL_NAME>().getTrimmedString(0),
|
||||
handlerContext.currentStep);
|
||||
|
||||
if (well_names.empty())
|
||||
this->invalidNamePattern(record.getItem<WELL_NAME>().getTrimmedString(0),
|
||||
handlerContext.currentStep, parseContext,
|
||||
errors, handlerContext.keyword);
|
||||
|
||||
const auto rawProdIndex = record.getItem<PI>().get<double>(0);
|
||||
for (const auto& well_name : well_names) {
|
||||
// All wells in a single record *hopefully* have the same preferred phase...
|
||||
const auto& well = this->getWell(well_name, handlerContext.currentStep);
|
||||
const auto unitPI = (well.getPreferredPhase() == Phase::GAS) ? gasPI : liqPI;
|
||||
|
||||
auto well2 = std::make_shared<Well>(well);
|
||||
if (well2->updateWellProductivityIndex(usys.to_si(unitPI, rawProdIndex)))
|
||||
this->updateWell(std::move(well2), handlerContext.currentStep);
|
||||
|
||||
this->addWellGroupEvent(well_name, ScheduleEvents::WELL_PRODUCTIVITY_INDEX, handlerContext.currentStep);
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
void Schedule::handleWELSEGS(const HandlerContext& handlerContext, const ParseContext&, ErrorGuard&) {
|
||||
const auto& record1 = handlerContext.keyword.getRecord(0);
|
||||
const auto& wname = record1.getItem("WELL").getTrimmedString(0);
|
||||
@ -1735,6 +1778,7 @@ namespace {
|
||||
{ "WECON" , &Schedule::handleWECON },
|
||||
{ "WEFAC" , &Schedule::handleWEFAC },
|
||||
{ "WELOPEN" , &Schedule::handleWELOPEN },
|
||||
{ "WELPI" , &Schedule::handleWELPI },
|
||||
{ "WELSEGS" , &Schedule::handleWELSEGS },
|
||||
{ "WELSPECS", &Schedule::handleWELSPECS },
|
||||
{ "WELTARG" , &Schedule::handleWELTARG },
|
||||
|
@ -488,6 +488,15 @@ bool Well::updateInjection(std::shared_ptr<WellInjectionProperties> injection_ar
|
||||
return false;
|
||||
}
|
||||
|
||||
bool Well::updateWellProductivityIndex(const double prodIndex) {
|
||||
const auto update = this->productivity_index != prodIndex;
|
||||
if (update)
|
||||
this->productivity_index = prodIndex;
|
||||
|
||||
// Note order here: We must always run prepareWellPIScaling(), but that operation
|
||||
// *may* not lead to requiring updating the well state, so return 'update' if not.
|
||||
return this->connections->prepareWellPIScaling() || update;
|
||||
}
|
||||
|
||||
bool Well::updateHasProduced() {
|
||||
if (this->wtype.producer() && this->status == Status::OPEN) {
|
||||
@ -801,6 +810,21 @@ void Well::setInsertIndex(std::size_t index) {
|
||||
this->insert_index = index;
|
||||
}
|
||||
|
||||
void Well::applyWellProdIndexScaling(const double currentEffectivePI) {
|
||||
if (this->connections->empty())
|
||||
// No connections for this well. Unexpected.
|
||||
return;
|
||||
|
||||
if (this->productivity_index < 0.0)
|
||||
// WELPI not activated. Nothing to do.
|
||||
return;
|
||||
|
||||
if (this->productivity_index == currentEffectivePI)
|
||||
// No change in scaling.
|
||||
return;
|
||||
|
||||
this->connections->applyWellPIScaling(this->productivity_index / currentEffectivePI);
|
||||
}
|
||||
|
||||
const WellConnections& Well::getConnections() const {
|
||||
return *this->connections;
|
||||
@ -1486,6 +1510,10 @@ bool Well::operator==(const Well& data) const {
|
||||
this->getFoamProperties() == data.getFoamProperties() &&
|
||||
this->getStatus() == data.getStatus() &&
|
||||
this->guide_rate == data.guide_rate &&
|
||||
this->solvent_fraction == data.solvent_fraction &&
|
||||
this->hasProduced() == data.hasProduced() &&
|
||||
this->predictionMode() == data.predictionMode() &&
|
||||
this->productivity_index == data.productivity_index &&
|
||||
this->getTracerProperties() == data.getTracerProperties() &&
|
||||
this->getProductionProperties() == data.getProductionProperties() &&
|
||||
this->getInjectionProperties() == data.getInjectionProperties();
|
||||
|
@ -186,6 +186,36 @@ inline std::array< size_t, 3> directionIndices(const Opm::Connection::Direction
|
||||
return out;
|
||||
}
|
||||
|
||||
bool WellConnections::prepareWellPIScaling()
|
||||
{
|
||||
// New WELPI adjustment applies to all connections.
|
||||
auto update = ! (this->m_hasWellPIAdjustment && this->m_wellPIConnections.empty());
|
||||
|
||||
this->m_hasWellPIAdjustment = true;
|
||||
this->m_wellPIConnections.clear();
|
||||
|
||||
return update;
|
||||
}
|
||||
|
||||
void WellConnections::applyWellPIScaling(const double scaleFactor)
|
||||
{
|
||||
if (! this->m_hasWellPIAdjustment)
|
||||
return;
|
||||
|
||||
if (this->m_wellPIConnections.empty())
|
||||
// Apply to all connections
|
||||
for (auto& conn : this->m_connections)
|
||||
conn.scaleWellPi(scaleFactor);
|
||||
else {
|
||||
// Apply to active subset
|
||||
const auto nConn = std::min(this->m_wellPIConnections.size(),
|
||||
this->m_connections.size());
|
||||
|
||||
for (auto conn = 0*nConn; conn < nConn; ++conn)
|
||||
if (this->m_wellPIConnections[conn])
|
||||
this->m_connections[conn].scaleWellPi(scaleFactor);
|
||||
}
|
||||
}
|
||||
|
||||
void WellConnections::addConnection(int i, int j , int k ,
|
||||
std::size_t global_index,
|
||||
@ -413,6 +443,8 @@ inline std::array< size_t, 3> directionIndices(const Opm::Connection::Direction
|
||||
depth,
|
||||
css_ind,
|
||||
*perf_range);
|
||||
|
||||
this->excludeFromWellPI(std::distance(this->m_connections.begin(), prev));
|
||||
}
|
||||
}
|
||||
}
|
||||
@ -472,6 +504,7 @@ inline std::array< size_t, 3> directionIndices(const Opm::Connection::Direction
|
||||
|
||||
void WellConnections::add( Connection connection ) {
|
||||
m_connections.emplace_back( connection );
|
||||
this->excludeFromWellPI(this->m_connections.size() - 1);
|
||||
}
|
||||
|
||||
bool WellConnections::allConnectionsShut( ) const {
|
||||
@ -533,6 +566,32 @@ inline std::array< size_t, 3> directionIndices(const Opm::Connection::Direction
|
||||
}
|
||||
}
|
||||
|
||||
void WellConnections::excludeFromWellPI(const std::size_t connID) {
|
||||
assert (!this->m_connections.empty() &&
|
||||
"Connection set must be non-empty before calling excludeFromWellPI");
|
||||
|
||||
if (!this->m_hasWellPIAdjustment)
|
||||
// Connection set not prepared for WELPI. Common case. Nothing to do.
|
||||
return;
|
||||
|
||||
if (connID >= this->m_connections.size())
|
||||
throw std::invalid_argument {
|
||||
"Cannot exclude connection ID outside known range. "
|
||||
"Expected 0.." + std::to_string(this->m_connections.size() - 1)
|
||||
+ ", but got " + std::to_string(connID)
|
||||
};
|
||||
|
||||
if (this->m_wellPIConnections.empty())
|
||||
// WELPI applies to all connections. Prepare to exclude 'connID'.
|
||||
this->m_wellPIConnections.assign(this->m_connections.size(), true);
|
||||
|
||||
if (this->m_wellPIConnections.size() < this->m_connections.size())
|
||||
// WELPI applies to subset of connections. Must also exclude those
|
||||
// that are not already in set of explicitly included connections.
|
||||
this->m_wellPIConnections.resize(this->m_connections.size(), false);
|
||||
|
||||
this->m_wellPIConnections[connID] = false;
|
||||
}
|
||||
|
||||
size_t WellConnections::findClosestConnection(int oi, int oj, double oz, size_t start_pos)
|
||||
{
|
||||
@ -566,6 +625,8 @@ inline std::array< size_t, 3> directionIndices(const Opm::Connection::Direction
|
||||
bool WellConnections::operator==( const WellConnections& rhs ) const {
|
||||
return this->size() == rhs.size() &&
|
||||
this->m_ordering == rhs.m_ordering &&
|
||||
this->m_hasWellPIAdjustment == rhs.m_hasWellPIAdjustment &&
|
||||
this->m_wellPIConnections == rhs.m_wellPIConnections &&
|
||||
std::equal( this->begin(), this->end(), rhs.begin() );
|
||||
}
|
||||
|
||||
@ -573,12 +634,60 @@ inline std::array< size_t, 3> directionIndices(const Opm::Connection::Direction
|
||||
return !( *this == rhs );
|
||||
}
|
||||
|
||||
namespace {
|
||||
template <typename S1FwdIter, typename S2FwdIter, typename Predicate>
|
||||
std::pair<S1FwdIter, S2FwdIter>
|
||||
remove_if(S1FwdIter s1begin, S1FwdIter s1end, S2FwdIter s2begin, Predicate&& predicate)
|
||||
{
|
||||
auto ret = std::make_pair(s1begin, s2begin);
|
||||
|
||||
for (; s1begin != s1end; ++s1begin, ++s2begin) {
|
||||
if (predicate(*s1begin))
|
||||
// Skip those elements that match the predicate.
|
||||
continue;
|
||||
|
||||
if (ret.first != s1begin) {
|
||||
// Do nothing if src == dest, i.e. if we've not
|
||||
// skipped any sequence elements. Otherwise, move
|
||||
// down.
|
||||
*ret.first = std::move(*s1begin);
|
||||
*ret.second = *s2begin;
|
||||
}
|
||||
|
||||
++ret.first;
|
||||
++ret.second;
|
||||
}
|
||||
|
||||
return ret;
|
||||
}
|
||||
}
|
||||
|
||||
void WellConnections::filter(const ActiveGridCells& grid) {
|
||||
auto new_end = std::remove_if(m_connections.begin(),
|
||||
m_connections.end(),
|
||||
[&grid](const Connection& c) { return !grid.cellActive(c.getI(), c.getJ(), c.getK()); });
|
||||
m_connections.erase(new_end, m_connections.end());
|
||||
auto isInactive = [&grid](const Connection& c) {
|
||||
return !grid.cellActive(c.getI(), c.getJ(), c.getK());
|
||||
};
|
||||
|
||||
if (!this->m_hasWellPIAdjustment || this->m_wellPIConnections.empty()) {
|
||||
// Common case. Either no WELPI or WELPI applies to all connections.
|
||||
// Don't need to take WELPI scaling subset into account.
|
||||
auto new_end = std::remove_if(m_connections.begin(), m_connections.end(), isInactive);
|
||||
m_connections.erase(new_end, m_connections.end());
|
||||
}
|
||||
else {
|
||||
// Special case. Subset of connections subject to WELPI scaling.
|
||||
// Preserve those too.
|
||||
if (this->m_wellPIConnections.size() < this->m_connections.size())
|
||||
this->m_wellPIConnections.resize(this->m_connections.size(), false);
|
||||
|
||||
auto new_end = remove_if(this->m_connections.begin(), this->m_connections.end(),
|
||||
this->m_wellPIConnections.begin(), isInactive);
|
||||
|
||||
if (new_end.first != this->m_connections.end()) {
|
||||
// Prune moved from sequence elements.
|
||||
this->m_connections.erase(new_end.first, this->m_connections.end());
|
||||
this->m_wellPIConnections.erase(new_end.second, this->m_wellPIConnections.end());
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
|
@ -43,7 +43,15 @@
|
||||
#include <opm/parser/eclipse/EclipseState/Runspec.hpp>
|
||||
#include <opm/parser/eclipse/EclipseState/Schedule/Schedule.hpp>
|
||||
|
||||
#include <opm/parser/eclipse/Units/Units.hpp>
|
||||
|
||||
namespace {
|
||||
double cp_rm3_per_db()
|
||||
{
|
||||
return Opm::prefix::centi*Opm::unit::Poise * Opm::unit::cubic(Opm::unit::meter)
|
||||
/ (Opm::unit::day * Opm::unit::barsa);
|
||||
}
|
||||
|
||||
Opm::WellConnections loadCOMPDAT(const std::string& compdat_keyword) {
|
||||
Opm::EclipseGrid grid(10,10,10);
|
||||
Opm::TableManager tables;
|
||||
@ -340,3 +348,149 @@ BOOST_AUTO_TEST_CASE(loadCOMPDATTESTSPE9) {
|
||||
BOOST_CHECK_MESSAGE( !conn.ctfAssignedFromInput(), "Calculated SPE9 CTF values must NOT be assigned from input");
|
||||
}
|
||||
}
|
||||
|
||||
BOOST_AUTO_TEST_CASE(ApplyWellPI) {
|
||||
const auto deck = Opm::Parser{}.parseString(R"(RUNSPEC
|
||||
DIMENS
|
||||
10 10 3 /
|
||||
|
||||
START
|
||||
5 OCT 2020 /
|
||||
|
||||
GRID
|
||||
DXV
|
||||
10*100 /
|
||||
DYV
|
||||
10*100 /
|
||||
DZV
|
||||
3*10 /
|
||||
DEPTHZ
|
||||
121*2000 /
|
||||
|
||||
ACTNUM
|
||||
100*1
|
||||
99*1 0
|
||||
100*1
|
||||
/
|
||||
|
||||
PERMX
|
||||
300*100 /
|
||||
PERMY
|
||||
300*100 /
|
||||
PERMZ
|
||||
300*100 /
|
||||
PORO
|
||||
300*0.3 /
|
||||
|
||||
SCHEDULE
|
||||
WELSPECS
|
||||
'P' 'G' 10 10 2005 'LIQ' /
|
||||
/
|
||||
|
||||
COMPDAT
|
||||
'P' 0 0 1 3 OPEN 1 100 /
|
||||
/
|
||||
|
||||
TSTEP
|
||||
10
|
||||
/
|
||||
|
||||
END
|
||||
)");
|
||||
|
||||
const auto es = Opm::EclipseState{ deck };
|
||||
const auto sched = Opm::Schedule{ deck, es };
|
||||
|
||||
const auto expectCF = 100.0*cp_rm3_per_db();
|
||||
|
||||
auto connP = sched.getWell("P", 0).getConnections();
|
||||
for (const auto& conn : connP) {
|
||||
BOOST_CHECK_CLOSE(conn.CF(), expectCF, 1.0e-10);
|
||||
}
|
||||
|
||||
connP.applyWellPIScaling(2.0); // No "prepare" -> no change.
|
||||
for (const auto& conn : connP) {
|
||||
BOOST_CHECK_CLOSE(conn.CF(), expectCF, 1.0e-10);
|
||||
}
|
||||
|
||||
// All CFs scaled by factor 2.
|
||||
BOOST_CHECK_MESSAGE( connP.prepareWellPIScaling(), "First call to prepareWellPIScaling must be a state change");
|
||||
BOOST_CHECK_MESSAGE(!connP.prepareWellPIScaling(), "Second call to prepareWellPIScaling must NOT be a state change");
|
||||
connP.applyWellPIScaling(2.0);
|
||||
for (const auto& conn : connP) {
|
||||
BOOST_CHECK_CLOSE(conn.CF(), 2.0*expectCF, 1.0e-10);
|
||||
}
|
||||
|
||||
// Reset CF -- simulating COMPDAT record (inactive cell)
|
||||
connP.addConnection(9, 9, 1, // 10, 10, 2
|
||||
199,
|
||||
2015.0,
|
||||
Opm::Connection::State::OPEN,
|
||||
50.0*cp_rm3_per_db(),
|
||||
0.123,
|
||||
0.234,
|
||||
0.157,
|
||||
0.0,
|
||||
1);
|
||||
|
||||
BOOST_REQUIRE_EQUAL(connP.size(), std::size_t{3});
|
||||
|
||||
BOOST_CHECK_CLOSE(connP[0].CF(), 2.0*expectCF , 1.0e-10);
|
||||
BOOST_CHECK_CLOSE(connP[1].CF(), 2.0*expectCF , 1.0e-10);
|
||||
BOOST_CHECK_CLOSE(connP[2].CF(), 50.0*cp_rm3_per_db(), 1.0e-10);
|
||||
|
||||
// Should not apply to connection whose CF was manually specified
|
||||
connP.applyWellPIScaling(2.0);
|
||||
|
||||
BOOST_CHECK_CLOSE(connP[0].CF(), 4.0*expectCF , 1.0e-10);
|
||||
BOOST_CHECK_CLOSE(connP[1].CF(), 4.0*expectCF , 1.0e-10);
|
||||
BOOST_CHECK_CLOSE(connP[2].CF(), 50.0*cp_rm3_per_db(), 1.0e-10);
|
||||
|
||||
// Prepare new scaling. Simulating new WELPI record.
|
||||
// New scaling applies to all connections.
|
||||
BOOST_CHECK_MESSAGE(connP.prepareWellPIScaling(), "Third call to prepareWellPIScaling must be a state change");
|
||||
connP.applyWellPIScaling(2.0);
|
||||
|
||||
BOOST_CHECK_CLOSE(connP[0].CF(), 8.0*expectCF , 1.0e-10);
|
||||
BOOST_CHECK_CLOSE(connP[1].CF(), 8.0*expectCF , 1.0e-10);
|
||||
BOOST_CHECK_CLOSE(connP[2].CF(), 100.0*cp_rm3_per_db(), 1.0e-10);
|
||||
|
||||
// Reset CF -- simulating COMPDAT record (active cell)
|
||||
connP.addConnection(8, 9, 1, // 10, 10, 2
|
||||
198,
|
||||
2015.0,
|
||||
Opm::Connection::State::OPEN,
|
||||
50.0*cp_rm3_per_db(),
|
||||
0.123,
|
||||
0.234,
|
||||
0.157,
|
||||
0.0,
|
||||
1);
|
||||
|
||||
BOOST_REQUIRE_EQUAL(connP.size(), std::size_t{4});
|
||||
connP.applyWellPIScaling(2.0);
|
||||
|
||||
BOOST_CHECK_CLOSE(connP[0].CF(), 16.0*expectCF , 1.0e-10);
|
||||
BOOST_CHECK_CLOSE(connP[1].CF(), 16.0*expectCF , 1.0e-10);
|
||||
BOOST_CHECK_CLOSE(connP[2].CF(), 200.0*cp_rm3_per_db(), 1.0e-10);
|
||||
BOOST_CHECK_CLOSE(connP[3].CF(), 50.0*cp_rm3_per_db(), 1.0e-10);
|
||||
|
||||
const auto& grid = es.getInputGrid();
|
||||
const auto actCells = Opm::ActiveGridCells {
|
||||
std::size_t{10}, std::size_t{10}, std::size_t{3},
|
||||
grid.getActiveMap().data(),
|
||||
grid.getNumActive()
|
||||
};
|
||||
|
||||
connP.filter(actCells);
|
||||
|
||||
BOOST_REQUIRE_EQUAL(connP.size(), std::size_t{3});
|
||||
BOOST_CHECK_CLOSE(connP[0].CF(), 16.0*expectCF , 1.0e-10);
|
||||
BOOST_CHECK_CLOSE(connP[1].CF(), 16.0*expectCF , 1.0e-10);
|
||||
BOOST_CHECK_CLOSE(connP[2].CF(), 50.0*cp_rm3_per_db(), 1.0e-10);
|
||||
|
||||
connP.applyWellPIScaling(2.0);
|
||||
BOOST_CHECK_CLOSE(connP[0].CF(), 32.0*expectCF , 1.0e-10);
|
||||
BOOST_CHECK_CLOSE(connP[1].CF(), 32.0*expectCF , 1.0e-10);
|
||||
BOOST_CHECK_CLOSE(connP[2].CF(), 50.0*cp_rm3_per_db(), 1.0e-10);
|
||||
}
|
||||
|
@ -58,6 +58,18 @@
|
||||
|
||||
using namespace Opm;
|
||||
|
||||
namespace {
|
||||
double liquid_PI_unit()
|
||||
{
|
||||
return UnitSystem::newMETRIC().to_si(UnitSystem::measure::liquid_productivity_index, 1.0);
|
||||
}
|
||||
|
||||
double cp_rm3_per_db()
|
||||
{
|
||||
return prefix::centi*unit::Poise * unit::cubic(unit::meter)
|
||||
/ (unit::day * unit::barsa);
|
||||
}
|
||||
}
|
||||
|
||||
static Schedule make_schedule(const std::string& deck_string) {
|
||||
const auto& deck = Parser{}.parseString(deck_string);
|
||||
@ -3702,3 +3714,82 @@ WLIFTOPT
|
||||
BOOST_CHECK(w3.alloc_extra_gas());
|
||||
}
|
||||
|
||||
BOOST_AUTO_TEST_CASE(WellPI) {
|
||||
const auto deck = Parser{}.parseString(R"(RUNSPEC
|
||||
START
|
||||
7 OCT 2020 /
|
||||
|
||||
DIMENS
|
||||
10 10 3 /
|
||||
|
||||
GRID
|
||||
DXV
|
||||
10*100.0 /
|
||||
DYV
|
||||
10*100.0 /
|
||||
DZV
|
||||
3*10.0 /
|
||||
|
||||
DEPTHZ
|
||||
121*2000.0 /
|
||||
|
||||
PERMX
|
||||
300*100.0 /
|
||||
PERMY
|
||||
300*100.0 /
|
||||
PERMZ
|
||||
300*10.0 /
|
||||
PORO
|
||||
300*0.3 /
|
||||
|
||||
SCHEDULE
|
||||
WELSPECS
|
||||
'P' 'G' 10 10 2005 'LIQ' /
|
||||
/
|
||||
COMPDAT
|
||||
'P' 0 0 1 3 OPEN 1 100 /
|
||||
/
|
||||
|
||||
TSTEP
|
||||
10
|
||||
/
|
||||
|
||||
WELPI
|
||||
'P' 200.0 /
|
||||
/
|
||||
|
||||
TSTEP
|
||||
10
|
||||
/
|
||||
|
||||
END
|
||||
)");
|
||||
|
||||
const auto es = EclipseState{ deck };
|
||||
const auto sched = Schedule{ deck, es };
|
||||
|
||||
// Apply WELPI before seeing WELPI data
|
||||
{
|
||||
const auto expectCF = 100.0*cp_rm3_per_db();
|
||||
auto wellP = sched.getWell("P", 0);
|
||||
|
||||
wellP.applyWellProdIndexScaling(2.7182818);
|
||||
for (const auto& conn : wellP.getConnections()) {
|
||||
BOOST_CHECK_CLOSE(conn.CF(), expectCF, 1.0e-10);
|
||||
}
|
||||
}
|
||||
|
||||
// Apply WELPI after seeing WELPI data.
|
||||
{
|
||||
const auto expectCF = (200.0 / 100.0) * 100.0*cp_rm3_per_db();
|
||||
auto wellP = sched.getWell("P", 1);
|
||||
|
||||
wellP.applyWellProdIndexScaling(100.0*liquid_PI_unit());
|
||||
for (const auto& conn : wellP.getConnections()) {
|
||||
BOOST_CHECK_CLOSE(conn.CF(), expectCF, 1.0e-10);
|
||||
}
|
||||
}
|
||||
|
||||
BOOST_CHECK_MESSAGE(sched.hasWellGroupEvent("P", ScheduleEvents::WELL_PRODUCTIVITY_INDEX, 1),
|
||||
"Must have WELL_PRODUCTIVITY_INDEX event at report step 1");
|
||||
}
|
||||
|
@ -49,14 +49,12 @@
|
||||
|
||||
using namespace Opm;
|
||||
|
||||
|
||||
namespace Opm {
|
||||
inline std::ostream& operator<<( std::ostream& stream, const Connection& c ) {
|
||||
return stream << "(" << c.getI() << "," << c.getJ() << "," << c.getK() << ")";
|
||||
}
|
||||
inline std::ostream& operator<<( std::ostream& stream, const Well& well ) {
|
||||
return stream << "(" << well.name() << ")";
|
||||
}
|
||||
namespace {
|
||||
double cp_rm3_per_db()
|
||||
{
|
||||
return prefix::centi*unit::Poise * unit::cubic(unit::meter)
|
||||
/ (unit::day * unit::barsa);
|
||||
}
|
||||
}
|
||||
|
||||
BOOST_AUTO_TEST_CASE(WellCOMPDATtestTRACK) {
|
||||
@ -1159,3 +1157,93 @@ WCONINJE
|
||||
}
|
||||
|
||||
|
||||
BOOST_AUTO_TEST_CASE(WellPI) {
|
||||
const auto deck = Parser{}.parseString(R"(RUNSPEC
|
||||
START
|
||||
7 OCT 2020 /
|
||||
|
||||
DIMENS
|
||||
10 10 3 /
|
||||
|
||||
GRID
|
||||
DXV
|
||||
10*100.0 /
|
||||
DYV
|
||||
10*100.0 /
|
||||
DZV
|
||||
3*10.0 /
|
||||
|
||||
DEPTHZ
|
||||
121*2000.0 /
|
||||
|
||||
PERMX
|
||||
300*100.0 /
|
||||
PERMY
|
||||
300*100.0 /
|
||||
PERMZ
|
||||
300*10.0 /
|
||||
PORO
|
||||
300*0.3 /
|
||||
|
||||
SCHEDULE
|
||||
WELSPECS
|
||||
'P' 'G' 10 10 2005 'LIQ' /
|
||||
/
|
||||
COMPDAT
|
||||
'P' 0 0 1 3 OPEN 1 100 /
|
||||
/
|
||||
|
||||
END
|
||||
)");
|
||||
|
||||
const auto es = EclipseState{ deck };
|
||||
const auto sched = Schedule{ deck, es };
|
||||
|
||||
const auto expectCF = 100.0*cp_rm3_per_db();
|
||||
|
||||
auto wellP = sched.getWell("P", 0);
|
||||
|
||||
for (const auto& conn : wellP.getConnections()) {
|
||||
BOOST_CHECK_CLOSE(conn.CF(), expectCF, 1.0e-10);
|
||||
}
|
||||
|
||||
// Simulate applying WELPI before WELPI keyword. No effect.
|
||||
wellP.applyWellProdIndexScaling(2.7182818);
|
||||
for (const auto& conn : wellP.getConnections()) {
|
||||
BOOST_CHECK_CLOSE(conn.CF(), expectCF, 1.0e-10);
|
||||
}
|
||||
|
||||
// Simulate applying WELPI after seeing
|
||||
//
|
||||
// WELPI
|
||||
// P 2 /
|
||||
// /
|
||||
//
|
||||
// (ignoring units of measure)
|
||||
BOOST_CHECK_MESSAGE( wellP.updateWellProductivityIndex(2.0), "First call to updateWellProductivityIndex() must be a state change");
|
||||
BOOST_CHECK_MESSAGE(!wellP.updateWellProductivityIndex(2.0), "Second call to updateWellProductivityIndex() must NOT be a state change");
|
||||
|
||||
// Want PI=2, but actual/effective PI=1 => scale CF by 2.0/1.0.
|
||||
wellP.applyWellProdIndexScaling(1.0);
|
||||
for (const auto& conn : wellP.getConnections()) {
|
||||
BOOST_CHECK_CLOSE(conn.CF(), 2.0*expectCF, 1.0e-10);
|
||||
}
|
||||
|
||||
// Repeated application of WELPI multiplies scaling factors.
|
||||
wellP.applyWellProdIndexScaling(1.0);
|
||||
for (const auto& conn : wellP.getConnections()) {
|
||||
BOOST_CHECK_CLOSE(conn.CF(), 4.0*expectCF, 1.0e-10);
|
||||
}
|
||||
|
||||
// New WELPI record does not reset the scaling factors
|
||||
wellP.updateWellProductivityIndex(3.0);
|
||||
for (const auto& conn : wellP.getConnections()) {
|
||||
BOOST_CHECK_CLOSE(conn.CF(), 4.0*expectCF, 1.0e-10);
|
||||
}
|
||||
|
||||
// Effective PI=desired PI => no scaling change
|
||||
wellP.applyWellProdIndexScaling(3.0);
|
||||
for (const auto& conn : wellP.getConnections()) {
|
||||
BOOST_CHECK_CLOSE(conn.CF(), 4.0*expectCF, 1.0e-10);
|
||||
}
|
||||
}
|
||||
|
Loading…
Reference in New Issue
Block a user