mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
Add all present controls to wells, not just the active one.
Some restructuring to support more proper well handling, esp. group control and surface rate controls in general.
This commit is contained in:
parent
432b9d4473
commit
d7512bdeb6
@ -25,6 +25,7 @@
|
|||||||
#include <opm/core/utility/ErrorMacros.hpp>
|
#include <opm/core/utility/ErrorMacros.hpp>
|
||||||
#include <opm/core/utility/Units.hpp>
|
#include <opm/core/utility/Units.hpp>
|
||||||
#include <opm/core/WellCollection.hpp>
|
#include <opm/core/WellCollection.hpp>
|
||||||
|
#include <opm/core/fluid/blackoil/phaseUsageFromDeck.hpp>
|
||||||
|
|
||||||
#include <tr1/array>
|
#include <tr1/array>
|
||||||
#include <cmath>
|
#include <cmath>
|
||||||
@ -37,10 +38,10 @@ namespace
|
|||||||
struct WellData
|
struct WellData
|
||||||
{
|
{
|
||||||
WellType type;
|
WellType type;
|
||||||
WellControlType control;
|
// WellControlType control;
|
||||||
double target;
|
// double target;
|
||||||
double reference_bhp_depth;
|
double reference_bhp_depth;
|
||||||
SurfaceComponent injected_phase;
|
// Opm::InjectionSpecification::InjectorType injected_phase;
|
||||||
};
|
};
|
||||||
|
|
||||||
|
|
||||||
@ -50,14 +51,18 @@ namespace
|
|||||||
double well_index;
|
double well_index;
|
||||||
};
|
};
|
||||||
|
|
||||||
|
namespace ProductionControl
|
||||||
int prod_control_mode(const std::string& control)
|
|
||||||
{
|
{
|
||||||
const int num_prod_control_modes = 8;
|
enum Mode { ORAT, WRAT, GRAT,
|
||||||
|
LRAT, CRAT, RESV,
|
||||||
|
BHP, THP, GRUP };
|
||||||
|
Mode mode(const std::string& control)
|
||||||
|
{
|
||||||
|
const int num_prod_control_modes = 9;
|
||||||
static std::string prod_control_modes[num_prod_control_modes] =
|
static std::string prod_control_modes[num_prod_control_modes] =
|
||||||
{std::string("ORAT"), std::string("WRAT"), std::string("GRAT"),
|
{std::string("ORAT"), std::string("WRAT"), std::string("GRAT"),
|
||||||
std::string("LRAT"), std::string("RESV"), std::string("BHP"),
|
std::string("LRAT"), std::string("CRAT"), std::string("RESV"),
|
||||||
std::string("THP"), std::string("GRUP") };
|
std::string("BHP"), std::string("THP"), std::string("GRUP") };
|
||||||
int m = -1;
|
int m = -1;
|
||||||
for (int i=0; i<num_prod_control_modes; ++i) {
|
for (int i=0; i<num_prod_control_modes; ++i) {
|
||||||
if (control == prod_control_modes[i]) {
|
if (control == prod_control_modes[i]) {
|
||||||
@ -66,14 +71,19 @@ namespace
|
|||||||
}
|
}
|
||||||
}
|
}
|
||||||
if (m >= 0) {
|
if (m >= 0) {
|
||||||
return m;
|
return static_cast<Mode>(m);
|
||||||
} else {
|
} else {
|
||||||
THROW("Unknown well control mode = " << control << " in input file");
|
THROW("Unknown well control mode = " << control << " in input file");
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
} // namespace ProductionControl
|
||||||
|
|
||||||
|
|
||||||
int inje_control_mode(const std::string& control)
|
namespace InjectionControl
|
||||||
|
{
|
||||||
|
enum Mode { RATE, RESV, BHP,
|
||||||
|
THP, GRUP };
|
||||||
|
Mode mode(const std::string& control)
|
||||||
{
|
{
|
||||||
const int num_inje_control_modes = 5;
|
const int num_inje_control_modes = 5;
|
||||||
static std::string inje_control_modes[num_inje_control_modes] =
|
static std::string inje_control_modes[num_inje_control_modes] =
|
||||||
@ -88,11 +98,12 @@ namespace
|
|||||||
}
|
}
|
||||||
|
|
||||||
if (m >= 0) {
|
if (m >= 0) {
|
||||||
return m;
|
return static_cast<Mode>(m);
|
||||||
} else {
|
} else {
|
||||||
THROW("Unknown well control mode = " << control << " in input file");
|
THROW("Unknown well control mode = " << control << " in input file");
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
} // namespace InjectionControl
|
||||||
|
|
||||||
|
|
||||||
std::tr1::array<double, 3> getCubeDim(const UnstructuredGrid& grid, int cell)
|
std::tr1::array<double, 3> getCubeDim(const UnstructuredGrid& grid, int cell)
|
||||||
@ -202,6 +213,9 @@ namespace Opm
|
|||||||
THROW("Needed field is missing in file");
|
THROW("Needed field is missing in file");
|
||||||
}
|
}
|
||||||
|
|
||||||
|
// Obtain phase usage data.
|
||||||
|
PhaseUsage pu = phaseUsageFromDeck(deck);
|
||||||
|
|
||||||
// These data structures will be filled in this constructor,
|
// These data structures will be filled in this constructor,
|
||||||
// then used to initialize the Wells struct.
|
// then used to initialize the Wells struct.
|
||||||
std::vector<std::string> well_names;
|
std::vector<std::string> well_names;
|
||||||
@ -298,6 +312,7 @@ namespace Opm
|
|||||||
|
|
||||||
// Set up reference depths that were defaulted. Count perfs.
|
// Set up reference depths that were defaulted. Count perfs.
|
||||||
int num_perfs = 0;
|
int num_perfs = 0;
|
||||||
|
ASSERT(grid.dimensions == 3);
|
||||||
for (int w = 0; w < num_wells; ++w) {
|
for (int w = 0; w < num_wells; ++w) {
|
||||||
num_perfs += wellperf_data[w].size();
|
num_perfs += wellperf_data[w].size();
|
||||||
if (well_data[w].reference_bhp_depth == -1e100) {
|
if (well_data[w].reference_bhp_depth == -1e100) {
|
||||||
@ -312,14 +327,37 @@ namespace Opm
|
|||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
// Get WCONINJE data
|
// Create the well data structures.
|
||||||
|
w_ = create_wells(pu.num_phases, num_wells, num_perfs);
|
||||||
|
if (!w_) {
|
||||||
|
THROW("Failed creating Wells struct.");
|
||||||
|
}
|
||||||
|
|
||||||
|
// Add wells.
|
||||||
|
for (int w = 0; w < num_wells; ++w) {
|
||||||
|
const int w_num_perf = wellperf_data[w].size();
|
||||||
|
std::vector<int> perf_cells(w_num_perf);
|
||||||
|
std::vector<double> perf_prodind(w_num_perf);
|
||||||
|
for (int perf = 0; perf < w_num_perf; ++perf) {
|
||||||
|
perf_cells[0] = wellperf_data[w][perf].cell;
|
||||||
|
perf_prodind[0] = wellperf_data[w][perf].well_index;
|
||||||
|
}
|
||||||
|
const double* comp_frac = NULL;
|
||||||
|
// We initialize all wells with a null component fraction,
|
||||||
|
// and must (for injection wells) overwrite it later.
|
||||||
|
add_well(well_data[w].type, well_data[w].reference_bhp_depth, w_num_perf,
|
||||||
|
comp_frac, &perf_cells[0], &perf_prodind[0], w_);
|
||||||
|
}
|
||||||
|
|
||||||
|
// Get WCONINJE data, add injection controls to wells.
|
||||||
if (deck.hasField("WCONINJE")) {
|
if (deck.hasField("WCONINJE")) {
|
||||||
const WCONINJE& wconinjes = deck.getWCONINJE();
|
const WCONINJE& wconinjes = deck.getWCONINJE();
|
||||||
const int num_wconinjes = wconinjes.wconinje.size();
|
const int num_wconinjes = wconinjes.wconinje.size();
|
||||||
for (int kw = 0; kw < num_wconinjes; ++kw) {
|
for (int kw = 0; kw < num_wconinjes; ++kw) {
|
||||||
|
const WconinjeLine& wci_line = wconinjes.wconinje[kw];
|
||||||
// Extract well name, or the part of the well name that
|
// Extract well name, or the part of the well name that
|
||||||
// comes before the '*'.
|
// comes before the '*'.
|
||||||
std::string name = wconinjes.wconinje[kw].well_;
|
std::string name = wci_line.well_;
|
||||||
std::string::size_type len = name.find('*');
|
std::string::size_type len = name.find('*');
|
||||||
if (len != std::string::npos) {
|
if (len != std::string::npos) {
|
||||||
name = name.substr(0, len);
|
name = name.substr(0, len);
|
||||||
@ -328,45 +366,63 @@ namespace Opm
|
|||||||
for (int wix = 0; wix < num_wells; ++wix) {
|
for (int wix = 0; wix < num_wells; ++wix) {
|
||||||
if (well_names[wix].compare(0,len, name) == 0) { //equal
|
if (well_names[wix].compare(0,len, name) == 0) { //equal
|
||||||
well_found = true;
|
well_found = true;
|
||||||
well_data[wix].type = INJECTOR;
|
ASSERT(well_data[wix].type == w_->type[wix]);
|
||||||
int m = inje_control_mode(wconinjes.wconinje[kw].control_mode_);
|
if (well_data[wix].type != INJECTOR) {
|
||||||
switch(m) {
|
THROW("Found WCONINJE entry for a non-injector well: " << well_names[wix]);
|
||||||
case 0: // RATE
|
|
||||||
well_data[wix].control = RATE;
|
|
||||||
well_data[wix].target = wconinjes.wconinje[kw].surface_flow_max_rate_;
|
|
||||||
break;
|
|
||||||
case 1: // RESV
|
|
||||||
well_data[wix].control = RATE;
|
|
||||||
well_data[wix].target = wconinjes.wconinje[kw].fluid_volume_max_rate_;
|
|
||||||
break;
|
|
||||||
case 2: // BHP
|
|
||||||
well_data[wix].control = BHP;
|
|
||||||
well_data[wix].target = wconinjes.wconinje[kw].BHP_limit_;
|
|
||||||
break;
|
|
||||||
case 3: // THP
|
|
||||||
well_data[wix].control = BHP;
|
|
||||||
well_data[wix].target = wconinjes.wconinje[kw].THP_limit_;
|
|
||||||
break;
|
|
||||||
case 4:
|
|
||||||
break;
|
|
||||||
default:
|
|
||||||
THROW("Unknown well control mode; WCONIJE = "
|
|
||||||
<< wconinjes.wconinje[kw].control_mode_
|
|
||||||
<< " in input file");
|
|
||||||
}
|
}
|
||||||
if (wconinjes.wconinje[kw].injector_type_ == "WATER") {
|
|
||||||
well_data[wix].injected_phase = WATER;
|
// Add all controls that are present in well.
|
||||||
} else if (wconinjes.wconinje[kw].injector_type_ == "OIL") {
|
int control_pos[5] = { -1, -1, -1, -1, -1 };
|
||||||
well_data[wix].injected_phase = OIL;
|
if (wci_line.surface_flow_max_rate_ > 0.0) {
|
||||||
} else if (wconinjes.wconinje[kw].injector_type_ == "GAS") {
|
control_pos[InjectionControl::RATE] = w_->ctrls[wix]->num;
|
||||||
well_data[wix].injected_phase = GAS;
|
const double distr[3] = { 1.0, 1.0, 1.0 };
|
||||||
} else {
|
append_well_controls(SURFACE_RATE, wci_line.surface_flow_max_rate_,
|
||||||
THROW("Error in injector specification, found no known fluid type.");
|
distr, wix, w_);
|
||||||
}
|
}
|
||||||
|
if (wci_line.reservoir_flow_max_rate_ > 0.0) {
|
||||||
|
control_pos[InjectionControl::RESV] = w_->ctrls[wix]->num;
|
||||||
|
const double distr[3] = { 1.0, 1.0, 1.0 };
|
||||||
|
append_well_controls(RESERVOIR_RATE, wci_line.reservoir_flow_max_rate_,
|
||||||
|
distr, wix, w_);
|
||||||
|
}
|
||||||
|
if (wci_line.BHP_limit_ > 0.0) {
|
||||||
|
control_pos[InjectionControl::BHP] = w_->ctrls[wix]->num;
|
||||||
|
append_well_controls(BHP, wci_line.BHP_limit_,
|
||||||
|
NULL, wix, w_);
|
||||||
|
}
|
||||||
|
if (wci_line.THP_limit_ > 0.0) {
|
||||||
|
THROW("We cannot handle THP limit for well " << well_names[wix]);
|
||||||
|
}
|
||||||
|
InjectionControl::Mode mode = InjectionControl::mode(wci_line.control_mode_);
|
||||||
|
const int cpos = control_pos[mode];
|
||||||
|
if (cpos == -1 && mode != InjectionControl::GRUP) {
|
||||||
|
THROW("Control mode type " << mode << " not present in well " << well_names[wix]);
|
||||||
|
}
|
||||||
|
set_current_control(wix, cpos, w_);
|
||||||
|
|
||||||
|
// Set well component fraction.
|
||||||
|
double cf[3] = { 0.0, 0.0, 0.0 };
|
||||||
|
if (wci_line.injector_type_ == "WATER") {
|
||||||
|
if (!pu.phase_used[BlackoilPhases::Aqua]) {
|
||||||
|
THROW("Water phase not used, yet found water-injecting well.");
|
||||||
|
}
|
||||||
|
cf[pu.phase_pos[BlackoilPhases::Aqua]] = 1.0;
|
||||||
|
} else if (wci_line.injector_type_ == "OIL") {
|
||||||
|
if (!pu.phase_used[BlackoilPhases::Liquid]) {
|
||||||
|
THROW("Oil phase not used, yet found oil-injecting well.");
|
||||||
|
}
|
||||||
|
cf[pu.phase_pos[BlackoilPhases::Liquid]] = 1.0;
|
||||||
|
} else if (wci_line.injector_type_ == "GAS") {
|
||||||
|
if (!pu.phase_used[BlackoilPhases::Vapour]) {
|
||||||
|
THROW("Water phase not used, yet found water-injecting well.");
|
||||||
|
}
|
||||||
|
cf[pu.phase_pos[BlackoilPhases::Vapour]] = 1.0;
|
||||||
|
}
|
||||||
|
std::copy(cf, cf + pu.num_phases, w_->comp_frac + wix*pu.num_phases);
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
if (!well_found) {
|
if (!well_found) {
|
||||||
THROW("Undefined well name: " << wconinjes.wconinje[kw].well_
|
THROW("Undefined well name: " << wci_line.well_
|
||||||
<< " in WCONINJE");
|
<< " in WCONINJE");
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
@ -376,61 +432,91 @@ namespace Opm
|
|||||||
if (deck.hasField("WCONPROD")) {
|
if (deck.hasField("WCONPROD")) {
|
||||||
const WCONPROD& wconprods = deck.getWCONPROD();
|
const WCONPROD& wconprods = deck.getWCONPROD();
|
||||||
const int num_wconprods = wconprods.wconprod.size();
|
const int num_wconprods = wconprods.wconprod.size();
|
||||||
std::cout << "num_wconprods = " <<num_wconprods << std::endl;
|
|
||||||
for (int kw = 0; kw < num_wconprods; ++kw) {
|
for (int kw = 0; kw < num_wconprods; ++kw) {
|
||||||
std::string name = wconprods.wconprod[kw].well_;
|
const WconprodLine& wcp_line = wconprods.wconprod[kw];
|
||||||
|
std::string name = wcp_line.well_;
|
||||||
std::string::size_type len = name.find('*');
|
std::string::size_type len = name.find('*');
|
||||||
if (len != std::string::npos) {
|
if (len != std::string::npos) {
|
||||||
name = name.substr(0, len);
|
name = name.substr(0, len);
|
||||||
}
|
}
|
||||||
|
|
||||||
bool well_found = false;
|
bool well_found = false;
|
||||||
for (int wix = 0; wix < num_wells; ++wix) {
|
for (int wix = 0; wix < num_wells; ++wix) {
|
||||||
if (well_names[wix].compare(0,len, name) == 0) { //equal
|
if (well_names[wix].compare(0,len, name) == 0) { //equal
|
||||||
well_found = true;
|
well_found = true;
|
||||||
well_data[wix].type = PRODUCER;
|
ASSERT(well_data[wix].type == w_->type[wix]);
|
||||||
int m = prod_control_mode(wconprods.wconprod[kw].control_mode_);
|
if (well_data[wix].type != PRODUCER) {
|
||||||
switch(m) {
|
THROW("Found WCONPROD entry for a non-producer well: " << well_names[wix]);
|
||||||
case 0: // ORAT
|
|
||||||
well_data[wix].control = RATE;
|
|
||||||
well_data[wix].target = wconprods.wconprod[kw].oil_max_rate_;
|
|
||||||
break;
|
|
||||||
case 1: // WRAT
|
|
||||||
well_data[wix].control = RATE;
|
|
||||||
well_data[wix].target = wconprods.wconprod[kw].water_max_rate_;
|
|
||||||
break;
|
|
||||||
case 2: // GRAT
|
|
||||||
well_data[wix].control = RATE;
|
|
||||||
well_data[wix].target = wconprods.wconprod[kw].gas_max_rate_;
|
|
||||||
break;
|
|
||||||
case 3: // LRAT
|
|
||||||
well_data[wix].control = RATE;
|
|
||||||
well_data[wix].target = wconprods.wconprod[kw].liquid_max_rate_;
|
|
||||||
break;
|
|
||||||
case 4: // RESV
|
|
||||||
well_data[wix].control = RATE;
|
|
||||||
well_data[wix].target = wconprods.wconprod[kw].fluid_volume_max_rate_;
|
|
||||||
break;
|
|
||||||
case 5: // BHP
|
|
||||||
well_data[wix].control = BHP;
|
|
||||||
well_data[wix].target = wconprods.wconprod[kw].BHP_limit_;
|
|
||||||
break;
|
|
||||||
case 6: // THP
|
|
||||||
well_data[wix].control = BHP;
|
|
||||||
well_data[wix].target = wconprods.wconprod[kw].THP_limit_;
|
|
||||||
break;
|
|
||||||
case 7:
|
|
||||||
// Handle group here.
|
|
||||||
break;
|
|
||||||
default:
|
|
||||||
THROW("Unknown well control mode; WCONPROD = "
|
|
||||||
<< wconprods.wconprod[kw].control_mode_
|
|
||||||
<< " in input file");
|
|
||||||
}
|
}
|
||||||
|
// Add all controls that are present in well.
|
||||||
|
int control_pos[9] = { -1, -1, -1, -1, -1, -1, -1, -1, -1 };
|
||||||
|
if (wcp_line.oil_max_rate_ > 0.0) {
|
||||||
|
if (!pu.phase_used[BlackoilPhases::Liquid]) {
|
||||||
|
THROW("Oil phase not active and ORAT control specified.");
|
||||||
|
}
|
||||||
|
control_pos[ProductionControl::ORAT] = w_->ctrls[wix]->num;
|
||||||
|
double distr[3] = { 0.0, 0.0, 0.0 };
|
||||||
|
distr[pu.phase_pos[BlackoilPhases::Liquid]] = 1.0;
|
||||||
|
append_well_controls(SURFACE_RATE, wcp_line.oil_max_rate_,
|
||||||
|
distr, wix, w_);
|
||||||
|
}
|
||||||
|
if (wcp_line.water_max_rate_ > 0.0) {
|
||||||
|
if (!pu.phase_used[BlackoilPhases::Aqua]) {
|
||||||
|
THROW("Water phase not active and WRAT control specified.");
|
||||||
|
}
|
||||||
|
control_pos[ProductionControl::WRAT] = w_->ctrls[wix]->num;
|
||||||
|
double distr[3] = { 0.0, 0.0, 0.0 };
|
||||||
|
distr[pu.phase_pos[BlackoilPhases::Aqua]] = 1.0;
|
||||||
|
append_well_controls(SURFACE_RATE, wcp_line.water_max_rate_,
|
||||||
|
distr, wix, w_);
|
||||||
|
}
|
||||||
|
if (wcp_line.gas_max_rate_ > 0.0) {
|
||||||
|
if (!pu.phase_used[BlackoilPhases::Vapour]) {
|
||||||
|
THROW("Gas phase not active and GRAT control specified.");
|
||||||
|
}
|
||||||
|
control_pos[ProductionControl::GRAT] = w_->ctrls[wix]->num;
|
||||||
|
double distr[3] = { 0.0, 0.0, 0.0 };
|
||||||
|
distr[pu.phase_pos[BlackoilPhases::Vapour]] = 1.0;
|
||||||
|
append_well_controls(SURFACE_RATE, wcp_line.gas_max_rate_,
|
||||||
|
distr, wix, w_);
|
||||||
|
}
|
||||||
|
if (wcp_line.liquid_max_rate_ > 0.0) {
|
||||||
|
if (!pu.phase_used[BlackoilPhases::Aqua]) {
|
||||||
|
THROW("Water phase not active and LRAT control specified.");
|
||||||
|
}
|
||||||
|
if (!pu.phase_used[BlackoilPhases::Liquid]) {
|
||||||
|
THROW("Oil phase not active and LRAT control specified.");
|
||||||
|
}
|
||||||
|
control_pos[ProductionControl::LRAT] = w_->ctrls[wix]->num;
|
||||||
|
double distr[3] = { 0.0, 0.0, 0.0 };
|
||||||
|
distr[pu.phase_pos[BlackoilPhases::Aqua]] = 1.0;
|
||||||
|
distr[pu.phase_pos[BlackoilPhases::Liquid]] = 1.0;
|
||||||
|
append_well_controls(SURFACE_RATE, wcp_line.liquid_max_rate_,
|
||||||
|
distr, wix, w_);
|
||||||
|
}
|
||||||
|
if (wcp_line.reservoir_flow_max_rate_ > 0.0) {
|
||||||
|
control_pos[ProductionControl::RESV] = w_->ctrls[wix]->num;
|
||||||
|
double distr[3] = { 1.0, 1.0, 1.0 };
|
||||||
|
append_well_controls(RESERVOIR_RATE, wcp_line.reservoir_flow_max_rate_,
|
||||||
|
distr, wix, w_);
|
||||||
|
}
|
||||||
|
if (wcp_line.BHP_limit_ > 0.0) {
|
||||||
|
control_pos[ProductionControl::BHP] = w_->ctrls[wix]->num;
|
||||||
|
append_well_controls(BHP, wcp_line.BHP_limit_,
|
||||||
|
NULL, wix, w_);
|
||||||
|
}
|
||||||
|
if (wcp_line.THP_limit_ > 0.0) {
|
||||||
|
THROW("We cannot handle THP limit for well " << well_names[wix]);
|
||||||
|
}
|
||||||
|
ProductionControl::Mode mode = ProductionControl::mode(wcp_line.control_mode_);
|
||||||
|
const int cpos = control_pos[mode];
|
||||||
|
if (cpos == -1 && mode != ProductionControl::GRUP) {
|
||||||
|
THROW("Control mode type " << mode << " not present in well " << well_names[wix]);
|
||||||
|
}
|
||||||
|
set_current_control(wix, cpos, w_);
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
if (!well_found) {
|
if (!well_found) {
|
||||||
THROW("Undefined well name: " << wconprods.wconprod[kw].well_
|
THROW("Undefined well name: " << wcp_line.well_
|
||||||
<< " in WCONPROD");
|
<< " in WCONPROD");
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
@ -438,6 +524,8 @@ namespace Opm
|
|||||||
|
|
||||||
// Get WELTARG data
|
// Get WELTARG data
|
||||||
if (deck.hasField("WELTARG")) {
|
if (deck.hasField("WELTARG")) {
|
||||||
|
THROW("We currently do not handle WELTARG.");
|
||||||
|
/*
|
||||||
const WELTARG& weltargs = deck.getWELTARG();
|
const WELTARG& weltargs = deck.getWELTARG();
|
||||||
const int num_weltargs = weltargs.weltarg.size();
|
const int num_weltargs = weltargs.weltarg.size();
|
||||||
for (int kw = 0; kw < num_weltargs; ++kw) {
|
for (int kw = 0; kw < num_weltargs; ++kw) {
|
||||||
@ -459,11 +547,13 @@ namespace Opm
|
|||||||
<< " in WELTARG");
|
<< " in WELTARG");
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
*/
|
||||||
}
|
}
|
||||||
|
|
||||||
// Debug output.
|
// Debug output.
|
||||||
#define EXTRA_OUTPUT
|
#define EXTRA_OUTPUT
|
||||||
#ifdef EXTRA_OUTPUT
|
#ifdef EXTRA_OUTPUT
|
||||||
|
/*
|
||||||
std::cout << "\t WELL DATA" << std::endl;
|
std::cout << "\t WELL DATA" << std::endl;
|
||||||
for(int i = 0; i< num_wells; ++i) {
|
for(int i = 0; i< num_wells; ++i) {
|
||||||
std::cout << i << ": " << well_data[i].type << " "
|
std::cout << i << ": " << well_data[i].type << " "
|
||||||
@ -478,19 +568,19 @@ namespace Opm
|
|||||||
<< wellperf_data[i][j].well_index << std::endl;
|
<< wellperf_data[i][j].well_index << std::endl;
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
*/
|
||||||
#endif
|
#endif
|
||||||
|
|
||||||
|
|
||||||
|
// Build the well_collection_ well group hierarchy.
|
||||||
if (deck.hasField("GRUPTREE")) {
|
if (deck.hasField("GRUPTREE")) {
|
||||||
std::cout << "Found gruptree" << std::endl;
|
std::cout << "Found gruptree" << std::endl;
|
||||||
const GRUPTREE& gruptree = deck.getGRUPTREE();
|
const GRUPTREE& gruptree = deck.getGRUPTREE();
|
||||||
|
|
||||||
std::map<std::string, std::string>::const_iterator it = gruptree.tree.begin();
|
std::map<std::string, std::string>::const_iterator it = gruptree.tree.begin();
|
||||||
for( ; it != gruptree.tree.end(); ++it) {
|
for( ; it != gruptree.tree.end(); ++it) {
|
||||||
well_collection_.addChild(it->first, it->second, deck);
|
well_collection_.addChild(it->first, it->second, deck);
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
for (size_t i = 0; i < welspecs.welspecs.size(); ++i) {
|
for (size_t i = 0; i < welspecs.welspecs.size(); ++i) {
|
||||||
WelspecsLine line = welspecs.welspecs[i];
|
WelspecsLine line = welspecs.welspecs[i];
|
||||||
well_collection_.addChild(line.name_, line.group_, deck);
|
well_collection_.addChild(line.name_, line.group_, deck);
|
||||||
@ -506,110 +596,92 @@ namespace Opm
|
|||||||
std::cout << well_collection_.getLeafNodes().size() << std::endl;
|
std::cout << well_collection_.getLeafNodes().size() << std::endl;
|
||||||
for (size_t i = 0; i < lines.size(); i++) {
|
for (size_t i = 0; i < lines.size(); i++) {
|
||||||
std::string name = lines[i].well_;
|
std::string name = lines[i].well_;
|
||||||
int index = well_names_to_index[name];
|
const int wix = well_names_to_index[name];
|
||||||
ASSERT(well_collection_.getLeafNodes()[index]->name() == name);
|
WellNode& wellnode = *well_collection_.getLeafNodes()[wix];
|
||||||
well_collection_.getLeafNodes()[index]->prodSpec().guide_rate_ = lines[i].guide_rate_;
|
ASSERT(wellnode.name() == name);
|
||||||
well_collection_.getLeafNodes()[index]->prodSpec().guide_rate_type_
|
if (well_data[wix].type == PRODUCER) {
|
||||||
= lines[i].phase_ == "OIL" ? ProductionSpecification::OIL : ProductionSpecification::RAT;
|
wellnode.prodSpec().guide_rate_ = lines[i].guide_rate_;
|
||||||
|
if (lines[i].phase_ == "OIL") {
|
||||||
|
wellnode.prodSpec().guide_rate_type_ = ProductionSpecification::OIL;
|
||||||
|
} else {
|
||||||
|
THROW("Guide rate type " << lines[i].phase_ << " specified for producer "
|
||||||
|
<< name << " in WGRUPCON, cannot handle.");
|
||||||
|
}
|
||||||
|
} else if (well_data[wix].type == INJECTOR) {
|
||||||
|
wellnode.injSpec().guide_rate_ = lines[i].guide_rate_;
|
||||||
|
if (lines[i].phase_ == "RAT") {
|
||||||
|
wellnode.injSpec().guide_rate_type_ = InjectionSpecification::RAT;
|
||||||
|
} else {
|
||||||
|
THROW("Guide rate type " << lines[i].phase_ << " specified for injector "
|
||||||
|
<< name << " in WGRUPCON, cannot handle.");
|
||||||
|
}
|
||||||
|
} else {
|
||||||
|
THROW("Unknown well type " << well_data[wix].type << " for well " << name);
|
||||||
|
}
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
well_collection_.calculateGuideRates();
|
well_collection_.calculateGuideRates();
|
||||||
|
|
||||||
// Apply guide rates:
|
// Apply guide rates:
|
||||||
for (size_t i = 0; i < well_data.size(); i++) {
|
for (size_t wix = 0; wix < well_data.size(); wix++) {
|
||||||
if (well_data[i].type == PRODUCER && (well_collection_.getLeafNodes()[i]->prodSpec().control_mode_ == ProductionSpecification::GRUP)) {
|
const WellNode& wellnode = *well_collection_.getLeafNodes()[wix];
|
||||||
switch (well_collection_.getLeafNodes()[i]->prodSpec().guide_rate_type_ ) {
|
if (well_data[wix].type == PRODUCER && (wellnode.prodSpec().control_mode_ == ProductionSpecification::GRUP)) {
|
||||||
|
ASSERT(w_->ctrls[wix]->current == -1);
|
||||||
|
switch (wellnode.prodSpec().guide_rate_type_ ) {
|
||||||
case ProductionSpecification::OIL:
|
case ProductionSpecification::OIL:
|
||||||
{
|
{
|
||||||
const ProductionSpecification& parent_prod_spec =
|
const ProductionSpecification& parent_prod_spec =
|
||||||
well_collection_.getLeafNodes()[i]->getParent()->prodSpec();
|
wellnode.getParent()->prodSpec();
|
||||||
double guide_rate = well_collection_.getLeafNodes()[i]->prodSpec().guide_rate_;
|
const double guide_rate = wellnode.prodSpec().guide_rate_;
|
||||||
well_data[i].target = guide_rate*parent_prod_spec.oil_max_rate_;
|
const double oil_target = guide_rate * parent_prod_spec.oil_max_rate_;
|
||||||
well_data[i].control = RATE;
|
double distr[3] = { 0.0, 0.0, 0.0 };
|
||||||
well_data[i].type = PRODUCER;
|
distr[pu.phase_pos[BlackoilPhases::Liquid]] = 1.0;
|
||||||
std::cout << "WARNING: Converting oil control to rate control!" << std::endl;
|
const int control_index = w_->ctrls[wix]->num;
|
||||||
break;
|
append_well_controls(SURFACE_RATE, oil_target, distr, wix, w_);
|
||||||
|
set_current_control(wix, control_index, w_);
|
||||||
}
|
}
|
||||||
case ProductionSpecification::NONE_GRT:
|
case ProductionSpecification::NONE_GRT:
|
||||||
{
|
{
|
||||||
// Will use the group control type:
|
// Will use the group control type:
|
||||||
const ProductionSpecification& parent_prod_spec =
|
const ProductionSpecification& parent_prod_spec =
|
||||||
well_collection_.getLeafNodes()[i]->getParent()->prodSpec();
|
wellnode.getParent()->prodSpec();
|
||||||
double guide_rate = well_collection_.getLeafNodes()[i]->prodSpec().guide_rate_;
|
const double guide_rate = wellnode.prodSpec().guide_rate_;
|
||||||
switch(parent_prod_spec.control_mode_) {
|
switch(parent_prod_spec.control_mode_) {
|
||||||
case ProductionSpecification::LRAT:
|
case ProductionSpecification::LRAT:
|
||||||
well_data[i].target = guide_rate * parent_prod_spec.liquid_max_rate_;
|
{
|
||||||
well_data[i].control = RATE;
|
const double liq_target = guide_rate * parent_prod_spec.liquid_max_rate_;
|
||||||
|
double distr[3] = { 0.0, 0.0, 0.0 };
|
||||||
|
distr[pu.phase_pos[BlackoilPhases::Aqua]] = 1.0;
|
||||||
|
distr[pu.phase_pos[BlackoilPhases::Liquid]] = 1.0;
|
||||||
|
append_well_controls(SURFACE_RATE, liq_target, distr, wix, w_);
|
||||||
break;
|
break;
|
||||||
|
}
|
||||||
default:
|
default:
|
||||||
THROW("Unhandled production specification control mode " << parent_prod_spec.control_mode_);
|
THROW("Unhandled parent production specification control mode " << parent_prod_spec.control_mode_);
|
||||||
break;
|
break;
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
default:
|
default:
|
||||||
THROW("Unhandled production specification guide rate type "
|
THROW("Unhandled production specification guide rate type "
|
||||||
<< well_collection_.getLeafNodes()[i]->prodSpec().guide_rate_type_);
|
<< wellnode.prodSpec().guide_rate_type_);
|
||||||
break;
|
break;
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
if (well_data[i].type == INJECTOR && (well_collection_.getLeafNodes()[i]->injSpec().control_mode_ == InjectionSpecification::GRUP)) {
|
if (well_data[wix].type == INJECTOR && (wellnode.injSpec().control_mode_ == InjectionSpecification::GRUP)) {
|
||||||
if (well_collection_.getLeafNodes()[i]->prodSpec().guide_rate_type_ == ProductionSpecification::RAT) {
|
ASSERT(w_->ctrls[wix]->current == -1);
|
||||||
well_data[i].injected_phase = WATER; // Default for now.
|
if (wellnode.injSpec().guide_rate_type_ == InjectionSpecification::RAT) {
|
||||||
well_data[i].control = RATE;
|
const double parent_surface_rate = wellnode.getParent()->injSpec().surface_flow_max_rate_;
|
||||||
well_data[i].type = INJECTOR;
|
const double guide_rate = wellnode.prodSpec().guide_rate_;
|
||||||
double parent_surface_rate = well_collection_.getLeafNodes()[i]->getParent()->injSpec().surface_flow_max_rate_;
|
const double target = guide_rate * parent_surface_rate;
|
||||||
double guide_rate = well_collection_.getLeafNodes()[i]->prodSpec().guide_rate_;
|
const double distr[3] = { 1.0, 1.0, 1.0 };
|
||||||
well_data[i].target = guide_rate * parent_surface_rate;
|
append_well_controls(SURFACE_RATE, target, distr, wix, w_);
|
||||||
std::cout << "Applying guide rate" << std::endl;
|
} else {
|
||||||
|
THROW("Unhandled injection specification guide rate type "
|
||||||
|
<< wellnode.injSpec().guide_rate_type_);
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
std::cout << "Making well structs" << std::endl;
|
|
||||||
// Set up the Wells struct.
|
|
||||||
w_ = create_wells(num_wells, num_perfs);
|
|
||||||
if (!w_) {
|
|
||||||
THROW("Failed creating Wells struct.");
|
|
||||||
}
|
|
||||||
double fracs[3][3] = { { 1.0, 0.0, 0.0 },
|
|
||||||
{ 0.0, 1.0, 0.0 },
|
|
||||||
{ 0.0, 0.0, 1.0 } };
|
|
||||||
for (int w = 0; w < num_wells; ++w) {
|
|
||||||
int nperf = wellperf_data[w].size();
|
|
||||||
std::vector<int> cells(nperf);
|
|
||||||
std::vector<double> wi(nperf);
|
|
||||||
for (int perf = 0; perf < nperf; ++perf) {
|
|
||||||
cells[perf] = wellperf_data[w][perf].cell;
|
|
||||||
wi[perf] = wellperf_data[w][perf].well_index;
|
|
||||||
}
|
|
||||||
const double* zfrac = (well_data[w].type == INJECTOR) ? fracs[well_data[w].injected_phase] : 0;
|
|
||||||
|
|
||||||
// DIRTY DIRTY HACK to temporarily make things work in spite of bugs in the deck reader.
|
|
||||||
if(well_data[w].type == INJECTOR && (well_data[w].injected_phase < 0 || well_data[w].injected_phase > 2)){
|
|
||||||
zfrac = fracs[WATER];
|
|
||||||
}
|
|
||||||
|
|
||||||
int ok = add_well(well_data[w].type, well_data[w].reference_bhp_depth, nperf,
|
|
||||||
zfrac, &cells[0], &wi[0], w_);
|
|
||||||
if (!ok) {
|
|
||||||
THROW("Failed to add a well.");
|
|
||||||
}
|
|
||||||
// We only append a single control at this point.
|
|
||||||
// TODO: Handle multiple controls.
|
|
||||||
if (well_data[w].type == PRODUCER && well_data[w].control == RATE) {
|
|
||||||
// Convention is that well rates for producers are negative.
|
|
||||||
well_data[w].target = -well_data[w].target;
|
|
||||||
}
|
|
||||||
ok = append_well_controls(well_data[w].control, well_data[w].target, w_->ctrls[w]);
|
|
||||||
w_->ctrls[w]->current = 0;
|
|
||||||
if (!ok) {
|
|
||||||
THROW("Failed to add well controls.");
|
|
||||||
}
|
|
||||||
}
|
|
||||||
|
|
||||||
std::cout << "Made well struct" << std::endl;
|
|
||||||
well_collection_.setWellsPointer(w_);
|
well_collection_.setWellsPointer(w_);
|
||||||
}
|
}
|
||||||
|
|
||||||
|
Loading…
Reference in New Issue
Block a user