mirror of
https://github.com/OPM/opm-simulators.git
synced 2024-12-26 01:01:00 -06:00
66b749c2eb
a few were probably forgotten.
848 lines
38 KiB
C++
848 lines
38 KiB
C++
/*
|
|
Copyright 2013 SINTEF ICT, Applied Mathematics.
|
|
Copyright 2014-2016 IRIS AS
|
|
Copyright 2015 Andreas Lauser
|
|
|
|
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/>.
|
|
*/
|
|
|
|
#include <utility>
|
|
#include <functional>
|
|
#include <algorithm>
|
|
#include <locale>
|
|
#include <opm/parser/eclipse/EclipseState/Schedule/Events.hpp>
|
|
#include <opm/core/utility/initHydroCarbonState.hpp>
|
|
#include <opm/core/well_controls.h>
|
|
#include <opm/core/wells/DynamicListEconLimited.hpp>
|
|
#include <opm/autodiff/BlackoilModel.hpp>
|
|
|
|
namespace Opm
|
|
{
|
|
|
|
template <class Implementation>
|
|
SimulatorBase<Implementation>::SimulatorBase(const ParameterGroup& param,
|
|
const Grid& grid,
|
|
DerivedGeology& geo,
|
|
BlackoilPropsAdFromDeck& props,
|
|
const RockCompressibility* rock_comp_props,
|
|
NewtonIterationBlackoilInterface& linsolver,
|
|
const double* gravity,
|
|
const bool has_disgas,
|
|
const bool has_vapoil,
|
|
std::shared_ptr<EclipseState> eclipse_state,
|
|
std::shared_ptr<Schedule> schedule,
|
|
std::shared_ptr<SummaryConfig> summary_config,
|
|
OutputWriter& output_writer,
|
|
const std::vector<double>& threshold_pressures_by_face,
|
|
const std::unordered_set<std::string>& defunct_well_names)
|
|
: param_(param),
|
|
model_param_(param),
|
|
solver_param_(param),
|
|
grid_(grid),
|
|
props_(props),
|
|
rock_comp_props_(rock_comp_props),
|
|
gravity_(gravity),
|
|
geo_(geo),
|
|
solver_(linsolver),
|
|
has_disgas_(has_disgas),
|
|
has_vapoil_(has_vapoil),
|
|
terminal_output_(param.getDefault("output_terminal", true)),
|
|
eclipse_state_(eclipse_state),
|
|
schedule_(schedule),
|
|
summary_config_(summary_config),
|
|
output_writer_(output_writer),
|
|
rateConverter_(props_.phaseUsage(), std::vector<int>(AutoDiffGrid::numCells(grid_), 0)),
|
|
threshold_pressures_by_face_(threshold_pressures_by_face),
|
|
is_parallel_run_( false ),
|
|
defunct_well_names_(defunct_well_names)
|
|
{
|
|
// Misc init.
|
|
const int num_cells = AutoDiffGrid::numCells(grid);
|
|
allcells_.resize(num_cells);
|
|
for (int cell = 0; cell < num_cells; ++cell) {
|
|
allcells_[cell] = cell;
|
|
}
|
|
#if HAVE_MPI
|
|
if ( solver_.parallelInformation().type() == typeid(ParallelISTLInformation) )
|
|
{
|
|
const ParallelISTLInformation& info =
|
|
boost::any_cast<const ParallelISTLInformation&>(solver_.parallelInformation());
|
|
// Only rank 0 does print to std::cout
|
|
terminal_output_ = terminal_output_ && ( info.communicator().rank() == 0 );
|
|
is_parallel_run_ = ( info.communicator().size() > 1 );
|
|
}
|
|
#endif
|
|
}
|
|
|
|
template <class Implementation>
|
|
SimulatorReport SimulatorBase<Implementation>::run(SimulatorTimer& timer,
|
|
ReservoirState& state)
|
|
{
|
|
WellState prev_well_state;
|
|
|
|
ExtraData extra;
|
|
if (output_writer_.isRestart()) {
|
|
// This is a restart, populate WellState and ReservoirState state objects from restart file
|
|
output_writer_.initFromRestartFile(props_.phaseUsage(), grid_, state, prev_well_state, extra);
|
|
initHydroCarbonState(state, props_.phaseUsage(), Opm::UgGridHelpers::numCells(grid_), has_disgas_, has_vapoil_);
|
|
initHysteresisParams(state);
|
|
}
|
|
|
|
// Create timers and file for writing timing info.
|
|
Opm::time::StopWatch solver_timer;
|
|
Opm::time::StopWatch step_timer;
|
|
Opm::time::StopWatch total_timer;
|
|
total_timer.start();
|
|
std::string tstep_filename = output_writer_.outputDirectory() + "/step_timing.txt";
|
|
std::ofstream tstep_os;
|
|
|
|
if ( output_writer_.output() ) {
|
|
if ( output_writer_.isIORank() )
|
|
tstep_os.open(tstep_filename.c_str());
|
|
}
|
|
|
|
// adaptive time stepping
|
|
const auto& events = schedule_->getEvents();
|
|
std::unique_ptr< AdaptiveTimeStepping > adaptiveTimeStepping;
|
|
if( param_.getDefault("timestep.adaptive", true ) )
|
|
{
|
|
|
|
if (param_.getDefault("use_TUNING", false)) {
|
|
adaptiveTimeStepping.reset( new AdaptiveTimeStepping( schedule_->getTuning(), timer.currentStepNum(), param_, terminal_output_ ) );
|
|
} else {
|
|
adaptiveTimeStepping.reset( new AdaptiveTimeStepping( param_, terminal_output_ ) );
|
|
}
|
|
if (output_writer_.isRestart()) {
|
|
if (extra.suggested_step > 0.0) {
|
|
adaptiveTimeStepping->setSuggestedNextStep(extra.suggested_step);
|
|
}
|
|
}
|
|
}
|
|
|
|
DynamicListEconLimited dynamic_list_econ_limited;
|
|
SimulatorReport report;
|
|
SimulatorReport stepReport;
|
|
|
|
bool ooip_computed = false;
|
|
std::vector<int> fipnum_global = eclipse_state_->get3DProperties().getIntGridProperty("FIPNUM").getData();
|
|
//Get compressed cell fipnum.
|
|
std::vector<int> fipnum(AutoDiffGrid::numCells(grid_));
|
|
if (fipnum_global.empty()) {
|
|
std::fill(fipnum.begin(), fipnum.end(), 0);
|
|
} else {
|
|
for (size_t c = 0; c < fipnum.size(); ++c) {
|
|
fipnum[c] = fipnum_global[AutoDiffGrid::globalCell(grid_)[c]];
|
|
}
|
|
}
|
|
std::vector<std::vector<double> > OOIP;
|
|
// Main simulation loop.
|
|
while (!timer.done()) {
|
|
// Report timestep.
|
|
step_timer.start();
|
|
if ( terminal_output_ )
|
|
{
|
|
std::ostringstream ss;
|
|
timer.report(ss);
|
|
OpmLog::note(ss.str());
|
|
}
|
|
|
|
// Create wells and well state.
|
|
WellsManager wells_manager(*eclipse_state_,
|
|
*schedule_,
|
|
timer.currentStepNum(),
|
|
Opm::UgGridHelpers::numCells(grid_),
|
|
Opm::UgGridHelpers::globalCell(grid_),
|
|
Opm::UgGridHelpers::cartDims(grid_),
|
|
Opm::UgGridHelpers::dimensions(grid_),
|
|
Opm::UgGridHelpers::cell2Faces(grid_),
|
|
Opm::UgGridHelpers::beginFaceCentroids(grid_),
|
|
dynamic_list_econ_limited,
|
|
is_parallel_run_,
|
|
defunct_well_names_);
|
|
const Wells* wells = wells_manager.c_wells();
|
|
WellState well_state;
|
|
well_state.initLegacy(wells, state, prev_well_state, props_.phaseUsage());
|
|
|
|
// give the polymer and surfactant simulators the chance to do their stuff
|
|
asImpl().handleAdditionalWellInflow(timer, wells_manager, well_state, wells);
|
|
|
|
// write the inital state at the report stage
|
|
if (timer.initialStep()) {
|
|
Dune::Timer perfTimer;
|
|
perfTimer.start();
|
|
|
|
// No per cell data is written for initial step, but will be
|
|
// for subsequent steps, when we have started simulating
|
|
output_writer_.writeTimeStepWithoutCellProperties( timer, state, well_state, {}, {} );
|
|
|
|
report.output_write_time += perfTimer.stop();
|
|
}
|
|
|
|
// Max oil saturation (for VPPARS), hysteresis update.
|
|
props_.updateSatOilMax(state.saturation());
|
|
props_.updateSatHyst(state.saturation(), allcells_);
|
|
|
|
// Compute reservoir volumes for RESV controls.
|
|
asImpl().computeRESV(timer.currentStepNum(), wells, state, well_state);
|
|
|
|
// Run a multiple steps of the solver depending on the time step control.
|
|
solver_timer.start();
|
|
|
|
const WellModel well_model(wells, &(wells_manager.wellCollection()), timer.currentStepNum());
|
|
|
|
std::unique_ptr<Solver> solver = asImpl().createSolver(well_model);
|
|
|
|
// Compute orignal FIP;
|
|
if (!ooip_computed) {
|
|
OOIP = solver->computeFluidInPlace(state, fipnum);
|
|
FIPUnitConvert(eclipse_state_->getUnits(), OOIP);
|
|
ooip_computed = true;
|
|
}
|
|
|
|
if( terminal_output_ )
|
|
{
|
|
std::ostringstream step_msg;
|
|
boost::posix_time::time_facet* facet = new boost::posix_time::time_facet("%d-%b-%Y");
|
|
step_msg.imbue(std::locale(std::locale::classic(), facet));
|
|
step_msg << "\nTime step " << std::setw(4) <<timer.currentStepNum()
|
|
<< " at day " << (double)unit::convert::to(timer.simulationTimeElapsed(), unit::day)
|
|
<< "/" << (double)unit::convert::to(timer.totalTime(), unit::day)
|
|
<< ", date = " << timer.currentDateTime();
|
|
OpmLog::info(step_msg.str());
|
|
}
|
|
|
|
// If sub stepping is enabled allow the solver to sub cycle
|
|
// in case the report steps are too large for the solver to converge
|
|
//
|
|
// \Note: The report steps are met in any case
|
|
// \Note: The sub stepping will require a copy of the state variables
|
|
if( adaptiveTimeStepping ) {
|
|
bool event = events.hasEvent(ScheduleEvents::NEW_WELL, timer.currentStepNum()) ||
|
|
events.hasEvent(ScheduleEvents::PRODUCTION_UPDATE, timer.currentStepNum()) ||
|
|
events.hasEvent(ScheduleEvents::INJECTION_UPDATE, timer.currentStepNum()) ||
|
|
events.hasEvent(ScheduleEvents::WELL_STATUS_CHANGE, timer.currentStepNum());
|
|
report += adaptiveTimeStepping->step( timer, *solver, state, well_state, event, output_writer_,
|
|
output_writer_.requireFIPNUM() ? &fipnum : nullptr );
|
|
}
|
|
else {
|
|
// solve for complete report step
|
|
stepReport = solver->step(timer, state, well_state);
|
|
report += stepReport;
|
|
|
|
if( terminal_output_ )
|
|
{
|
|
std::ostringstream iter_msg;
|
|
iter_msg << "Stepsize " << (double)unit::convert::to(timer.currentStepLength(), unit::day);
|
|
if (solver->wellIterations() != 0) {
|
|
iter_msg << " days well iterations = " << solver->wellIterations() << ", ";
|
|
}
|
|
iter_msg << "non-linear iterations = " << solver->nonlinearIterations()
|
|
<< ", total linear iterations = " << solver->linearIterations()
|
|
<< "\n";
|
|
OpmLog::info(iter_msg.str());
|
|
}
|
|
}
|
|
|
|
// update the derived geology (transmissibilities, pore volumes, etc) if the
|
|
// has geology changed for the next report step
|
|
const int nextTimeStepIdx = timer.currentStepNum() + 1;
|
|
if (nextTimeStepIdx < timer.numSteps()
|
|
&& events.hasEvent(ScheduleEvents::GEO_MODIFIER, nextTimeStepIdx)) {
|
|
// bring the contents of the keywords to the current state of the SCHEDULE
|
|
// section
|
|
//
|
|
// TODO (?): handle the parallel case (maybe this works out of the box)
|
|
const auto& miniDeck = schedule_->getModifierDeck(nextTimeStepIdx);
|
|
eclipse_state_->applyModifierDeck(miniDeck);
|
|
geo_.update(grid_, props_, *eclipse_state_, gravity_);
|
|
}
|
|
|
|
// take time that was used to solve system for this reportStep
|
|
solver_timer.stop();
|
|
|
|
// update timing.
|
|
report.solver_time += solver_timer.secsSinceStart();
|
|
|
|
// Compute current FIP.
|
|
std::vector<std::vector<double> > COIP;
|
|
COIP = solver->computeFluidInPlace(state, fipnum);
|
|
std::vector<double> OOIP_totals = FIPTotals(OOIP, state);
|
|
std::vector<double> COIP_totals = FIPTotals(COIP, state);
|
|
|
|
//Convert to correct units
|
|
FIPUnitConvert(eclipse_state_->getUnits(), COIP);
|
|
FIPUnitConvert(eclipse_state_->getUnits(), OOIP_totals);
|
|
FIPUnitConvert(eclipse_state_->getUnits(), COIP_totals);
|
|
|
|
if ( terminal_output_ )
|
|
{
|
|
outputFluidInPlace(OOIP_totals, COIP_totals,eclipse_state_->getUnits(), 0);
|
|
for (size_t reg = 0; reg < OOIP.size(); ++reg) {
|
|
outputFluidInPlace(OOIP[reg], COIP[reg], eclipse_state_->getUnits(), reg+1);
|
|
}
|
|
}
|
|
|
|
if ( terminal_output_ )
|
|
{
|
|
std::string msg;
|
|
msg = "Fully implicit solver took: " + std::to_string(stepReport.solver_time) + " seconds. Total solver time taken: " + std::to_string(report.solver_time) + " seconds.";
|
|
OpmLog::note(msg);
|
|
}
|
|
|
|
if ( tstep_os.is_open() ) {
|
|
stepReport.reportParam(tstep_os);
|
|
}
|
|
|
|
// Increment timer, remember well state.
|
|
++timer;
|
|
|
|
// write simulation state at the report stage
|
|
Dune::Timer perfTimer;
|
|
perfTimer.start();
|
|
const auto& physicalModel = solver->model();
|
|
output_writer_.writeTimeStep( timer, state, well_state, physicalModel );
|
|
report.output_write_time += perfTimer.stop();
|
|
|
|
prev_well_state = well_state;
|
|
|
|
asImpl().updateListEconLimited(solver, *schedule_, timer.currentStepNum(), wells,
|
|
well_state, dynamic_list_econ_limited);
|
|
}
|
|
|
|
// Stop timer and create timing report
|
|
total_timer.stop();
|
|
report.total_time = total_timer.secsSinceStart();
|
|
report.converged = true;
|
|
return report;
|
|
}
|
|
|
|
namespace SimFIBODetails {
|
|
typedef std::unordered_map<std::string, const Well* > WellMap;
|
|
|
|
inline WellMap
|
|
mapWells(const std::vector< const Well* >& wells)
|
|
{
|
|
WellMap wmap;
|
|
|
|
for (std::vector< const Well* >::const_iterator
|
|
w = wells.begin(), e = wells.end();
|
|
w != e; ++w)
|
|
{
|
|
wmap.insert(std::make_pair((*w)->name(), *w));
|
|
}
|
|
|
|
return wmap;
|
|
}
|
|
|
|
inline int
|
|
resv_control(const WellControls* ctrl)
|
|
{
|
|
int i, n = well_controls_get_num(ctrl);
|
|
|
|
bool match = false;
|
|
for (i = 0; (! match) && (i < n); ++i) {
|
|
match = well_controls_iget_type(ctrl, i) == RESERVOIR_RATE;
|
|
}
|
|
|
|
if (! match) { i = 0; }
|
|
|
|
return i - 1; // -1 if no match, undo final "++" otherwise
|
|
}
|
|
|
|
inline bool
|
|
is_resv(const Wells& wells,
|
|
const int w)
|
|
{
|
|
return (0 <= resv_control(wells.ctrls[w]));
|
|
}
|
|
|
|
inline bool
|
|
is_resv(const WellMap& wmap,
|
|
const std::string& name,
|
|
const std::size_t step)
|
|
{
|
|
bool match = false;
|
|
|
|
WellMap::const_iterator i = wmap.find(name);
|
|
|
|
if (i != wmap.end()) {
|
|
const Well* wp = i->second;
|
|
|
|
match = (wp->isProducer(step) &&
|
|
wp->getProductionProperties(step)
|
|
.hasProductionControl(WellProducer::RESV))
|
|
|| (wp->isInjector(step) &&
|
|
wp->getInjectionProperties(step)
|
|
.hasInjectionControl(WellInjector::RESV));
|
|
}
|
|
|
|
return match;
|
|
}
|
|
|
|
inline std::vector<int>
|
|
resvWells(const Wells* wells,
|
|
const std::size_t step,
|
|
const WellMap& wmap)
|
|
{
|
|
std::vector<int> resv_wells;
|
|
if( wells )
|
|
{
|
|
for (int w = 0, nw = wells->number_of_wells; w < nw; ++w) {
|
|
if (is_resv(*wells, w) ||
|
|
((wells->name[w] != 0) &&
|
|
is_resv(wmap, wells->name[w], step)))
|
|
{
|
|
resv_wells.push_back(w);
|
|
}
|
|
}
|
|
}
|
|
|
|
return resv_wells;
|
|
}
|
|
|
|
inline void
|
|
historyRates(const PhaseUsage& pu,
|
|
const WellProductionProperties& p,
|
|
std::vector<double>& rates)
|
|
{
|
|
assert (! p.predictionMode);
|
|
assert (rates.size() ==
|
|
std::vector<double>::size_type(pu.num_phases));
|
|
|
|
if (pu.phase_used[ BlackoilPhases::Aqua ]) {
|
|
const std::vector<double>::size_type
|
|
i = pu.phase_pos[ BlackoilPhases::Aqua ];
|
|
|
|
rates[i] = p.WaterRate;
|
|
}
|
|
|
|
if (pu.phase_used[ BlackoilPhases::Liquid ]) {
|
|
const std::vector<double>::size_type
|
|
i = pu.phase_pos[ BlackoilPhases::Liquid ];
|
|
|
|
rates[i] = p.OilRate;
|
|
}
|
|
|
|
if (pu.phase_used[ BlackoilPhases::Vapour ]) {
|
|
const std::vector<double>::size_type
|
|
i = pu.phase_pos[ BlackoilPhases::Vapour ];
|
|
|
|
rates[i] = p.GasRate;
|
|
}
|
|
}
|
|
} // namespace SimFIBODetails
|
|
|
|
template <class Implementation>
|
|
void SimulatorBase<Implementation>::handleAdditionalWellInflow(SimulatorTimer& /* timer */,
|
|
WellsManager& /* wells_manager */,
|
|
WellState& /* well_state */,
|
|
const Wells* /* wells */)
|
|
{ }
|
|
|
|
template <class Implementation>
|
|
auto SimulatorBase<Implementation>::createSolver(const WellModel& well_model)
|
|
-> std::unique_ptr<Solver>
|
|
{
|
|
auto model = std::unique_ptr<Model>(new Model(model_param_,
|
|
grid_,
|
|
props_,
|
|
geo_,
|
|
rock_comp_props_,
|
|
well_model,
|
|
solver_,
|
|
eclipse_state_,
|
|
schedule_,
|
|
summary_config_,
|
|
has_disgas_,
|
|
has_vapoil_,
|
|
terminal_output_));
|
|
|
|
if (!threshold_pressures_by_face_.empty()) {
|
|
model->setThresholdPressures(threshold_pressures_by_face_);
|
|
}
|
|
|
|
return std::unique_ptr<Solver>(new Solver(solver_param_, std::move(model)));
|
|
}
|
|
|
|
|
|
template <class Implementation>
|
|
void SimulatorBase<Implementation>::computeRESV(const std::size_t step,
|
|
const Wells* wells,
|
|
const BlackoilState& x,
|
|
WellState& xw)
|
|
{
|
|
typedef SimFIBODetails::WellMap WellMap;
|
|
|
|
const auto w_ecl = schedule_->getWells(step);
|
|
const WellMap& wmap = SimFIBODetails::mapWells(w_ecl);
|
|
|
|
const std::vector<int>& resv_wells = SimFIBODetails::resvWells(wells, step, wmap);
|
|
|
|
const std::size_t number_resv_wells = resv_wells.size();
|
|
std::size_t global_number_resv_wells = number_resv_wells;
|
|
#if HAVE_MPI
|
|
if ( solver_.parallelInformation().type() == typeid(ParallelISTLInformation) )
|
|
{
|
|
const auto& info =
|
|
boost::any_cast<const ParallelISTLInformation&>(solver_.parallelInformation());
|
|
global_number_resv_wells = info.communicator().sum(global_number_resv_wells);
|
|
if ( global_number_resv_wells )
|
|
{
|
|
// At least one process has resv wells. Therefore rate converter needs
|
|
// to calculate averages over regions that might cross process
|
|
// borders. This needs to be done by all processes and therefore
|
|
// outside of the next if statement.
|
|
rateConverter_.defineState(x, boost::any_cast<const ParallelISTLInformation&>(solver_.parallelInformation()));
|
|
}
|
|
}
|
|
else
|
|
#endif
|
|
{
|
|
if ( global_number_resv_wells )
|
|
{
|
|
rateConverter_.defineState(x);
|
|
}
|
|
}
|
|
|
|
if (! resv_wells.empty()) {
|
|
const PhaseUsage& pu = props_.phaseUsage();
|
|
const std::vector<double>::size_type np = props_.numPhases();
|
|
|
|
std::vector<double> distr (np);
|
|
std::vector<double> hrates(np);
|
|
std::vector<double> prates(np);
|
|
|
|
for (std::vector<int>::const_iterator
|
|
rp = resv_wells.begin(), e = resv_wells.end();
|
|
rp != e; ++rp)
|
|
{
|
|
WellControls* ctrl = wells->ctrls[*rp];
|
|
const bool is_producer = wells->type[*rp] == PRODUCER;
|
|
|
|
// RESV control mode, all wells
|
|
{
|
|
const int rctrl = SimFIBODetails::resv_control(ctrl);
|
|
|
|
if (0 <= rctrl) {
|
|
const std::vector<double>::size_type off = (*rp) * np;
|
|
|
|
if (is_producer) {
|
|
// Convert to positive rates to avoid issues
|
|
// in coefficient calculations.
|
|
std::transform(xw.wellRates().begin() + (off + 0*np),
|
|
xw.wellRates().begin() + (off + 1*np),
|
|
prates.begin(), std::negate<double>());
|
|
} else {
|
|
std::copy(xw.wellRates().begin() + (off + 0*np),
|
|
xw.wellRates().begin() + (off + 1*np),
|
|
prates.begin());
|
|
}
|
|
|
|
const int fipreg = 0; // Hack. Ignore FIP regions.
|
|
const int well_cell_top = wells->well_cells[wells->well_connpos[*rp]];
|
|
const int pvtreg = props_.cellPvtRegionIndex()[well_cell_top];
|
|
rateConverter_.calcCoeff(fipreg, pvtreg, distr);
|
|
|
|
well_controls_iset_distr(ctrl, rctrl, & distr[0]);
|
|
}
|
|
}
|
|
|
|
// RESV control, WCONHIST wells. A bit of duplicate
|
|
// work, regrettably.
|
|
if (is_producer && wells->name[*rp] != 0) {
|
|
WellMap::const_iterator i = wmap.find(wells->name[*rp]);
|
|
|
|
if (i != wmap.end()) {
|
|
const auto* wp = i->second;
|
|
|
|
const WellProductionProperties& p =
|
|
wp->getProductionProperties(step);
|
|
|
|
if (! p.predictionMode) {
|
|
// History matching (WCONHIST/RESV)
|
|
SimFIBODetails::historyRates(pu, p, hrates);
|
|
|
|
const int fipreg = 0; // Hack. Ignore FIP regions.
|
|
const int well_cell_top = wells->well_cells[wells->well_connpos[*rp]];
|
|
const int pvtreg = props_.cellPvtRegionIndex()[well_cell_top];
|
|
rateConverter_.calcCoeff(fipreg, pvtreg, distr);
|
|
|
|
// WCONHIST/RESV target is sum of all
|
|
// observed phase rates translated to
|
|
// reservoir conditions. Recall sign
|
|
// convention: Negative for producers.
|
|
const double target =
|
|
- std::inner_product(distr.begin(), distr.end(),
|
|
hrates.begin(), 0.0);
|
|
|
|
well_controls_clear(ctrl);
|
|
well_controls_assert_number_of_phases(ctrl, int(np));
|
|
|
|
static const double invalid_alq = -std::numeric_limits<double>::max();
|
|
static const int invalid_vfp = -std::numeric_limits<int>::max();
|
|
|
|
const int ok_resv =
|
|
well_controls_add_new(RESERVOIR_RATE, target,
|
|
invalid_alq, invalid_vfp,
|
|
& distr[0], ctrl);
|
|
|
|
// For WCONHIST the BHP limit is set to 1 atm.
|
|
// or a value specified using WELTARG
|
|
double bhp_limit = (p.BHPLimit > 0) ? p.BHPLimit : unit::convert::from(1.0, unit::atm);
|
|
const int ok_bhp =
|
|
well_controls_add_new(BHP, bhp_limit,
|
|
invalid_alq, invalid_vfp,
|
|
NULL, ctrl);
|
|
|
|
if (ok_resv != 0 && ok_bhp != 0) {
|
|
xw.currentControls()[*rp] = 0;
|
|
well_controls_set_current(ctrl, 0);
|
|
}
|
|
}
|
|
}
|
|
}
|
|
}
|
|
}
|
|
|
|
if( wells )
|
|
{
|
|
for (int w = 0, nw = wells->number_of_wells; w < nw; ++w) {
|
|
WellControls* ctrl = wells->ctrls[w];
|
|
const bool is_producer = wells->type[w] == PRODUCER;
|
|
if (!is_producer && wells->name[w] != 0) {
|
|
WellMap::const_iterator i = wmap.find(wells->name[w]);
|
|
if (i != wmap.end()) {
|
|
const auto* wp = i->second;
|
|
const WellInjectionProperties& injector = wp->getInjectionProperties(step);
|
|
if (!injector.predictionMode) {
|
|
//History matching WCONINJEH
|
|
static const double invalid_alq = -std::numeric_limits<double>::max();
|
|
static const int invalid_vfp = -std::numeric_limits<int>::max();
|
|
// For WCONINJEH the BHP limit is set to a large number
|
|
// or a value specified using WELTARG
|
|
double bhp_limit = (injector.BHPLimit > 0) ? injector.BHPLimit : std::numeric_limits<double>::max();
|
|
const int ok_bhp =
|
|
well_controls_add_new(BHP, bhp_limit,
|
|
invalid_alq, invalid_vfp,
|
|
NULL, ctrl);
|
|
if (!ok_bhp) {
|
|
OPM_THROW(std::runtime_error, "Failed to add well control.");
|
|
}
|
|
}
|
|
}
|
|
}
|
|
}
|
|
}
|
|
}
|
|
|
|
template <class Implementation>
|
|
void
|
|
SimulatorBase<Implementation>::FIPUnitConvert(const UnitSystem& units,
|
|
std::vector<std::vector<double> >& fip)
|
|
{
|
|
for (size_t i = 0; i < fip.size(); ++i) {
|
|
FIPUnitConvert(units, fip[i]);
|
|
}
|
|
}
|
|
|
|
template <class Implementation>
|
|
void
|
|
SimulatorBase<Implementation>::FIPUnitConvert(const UnitSystem& units, std::vector<double>& fip)
|
|
{
|
|
if (units.getType() == UnitSystem::UnitType::UNIT_TYPE_FIELD) {
|
|
fip[0] = unit::convert::to(fip[0], unit::stb);
|
|
fip[1] = unit::convert::to(fip[1], unit::stb);
|
|
fip[2] = unit::convert::to(fip[2], 1000*unit::cubic(unit::feet));
|
|
fip[3] = unit::convert::to(fip[3], 1000*unit::cubic(unit::feet));
|
|
fip[4] = unit::convert::to(fip[4], unit::stb);
|
|
fip[5] = unit::convert::to(fip[5], unit::stb);
|
|
fip[6] = unit::convert::to(fip[6], unit::psia);
|
|
}
|
|
else if (units.getType() == UnitSystem::UnitType::UNIT_TYPE_METRIC) {
|
|
fip[6] = unit::convert::to(fip[6], unit::barsa);
|
|
}
|
|
else {
|
|
OPM_THROW(std::runtime_error, "Unsupported unit type for fluid in place output.");
|
|
}
|
|
}
|
|
|
|
|
|
template <class Implementation>
|
|
std::vector<double>
|
|
SimulatorBase<Implementation>::FIPTotals(const std::vector<std::vector<double> >& fip, const ReservoirState& state)
|
|
{
|
|
std::vector<double> totals(7, 0.0);
|
|
for (int i = 0; i < 5; ++i) {
|
|
for (size_t reg = 0; reg < fip.size(); ++reg) {
|
|
totals[i] += fip[reg][i];
|
|
}
|
|
}
|
|
const int nc = Opm::AutoDiffGrid::numCells(grid_);
|
|
const int np = state.numPhases();
|
|
const PhaseUsage& pu = props_.phaseUsage();
|
|
const DataBlock s = Eigen::Map<const DataBlock>(& state.saturation()[0], nc, np);
|
|
std::vector<double> so(nc);
|
|
std::vector<double> sg(nc);
|
|
std::vector<double> hydrocarbon(nc);
|
|
// Using dummy indices if phase not used, the columns will not be accessed below if unused.
|
|
const int oilpos = pu.phase_used[BlackoilPhases::Liquid] ? pu.phase_pos[BlackoilPhases::Liquid] : 0;
|
|
const int gaspos = pu.phase_used[BlackoilPhases::Vapour] ? pu.phase_pos[BlackoilPhases::Vapour] : 0;
|
|
const auto& soCol = s.col(oilpos);
|
|
const auto& sgCol = s.col(gaspos);
|
|
for (unsigned c = 0; c < so.size(); ++ c) {
|
|
double mySo = 0.0;
|
|
if (pu.phase_used[BlackoilPhases::Liquid]) {
|
|
mySo = soCol[c];
|
|
}
|
|
double mySg = 0.0;
|
|
if (pu.phase_used[BlackoilPhases::Vapour]) {
|
|
mySg = sgCol[c];
|
|
}
|
|
so[c] = mySo;
|
|
sg[c] = mySg;
|
|
hydrocarbon[c] = mySo + mySg;
|
|
}
|
|
const std::vector<double> p = state.pressure();
|
|
if ( ! is_parallel_run_ )
|
|
{
|
|
double tmp = 0.0;
|
|
double tmp2 = 0.0;
|
|
for (unsigned i = 0; i < p.size(); ++i) {
|
|
tmp += p[i] * geo_.poreVolume()[i] * hydrocarbon[i];
|
|
tmp2 += geo_.poreVolume()[i] * hydrocarbon[i];
|
|
}
|
|
totals[5] = geo_.poreVolume().sum();
|
|
totals[6] = tmp/tmp2;
|
|
}
|
|
else
|
|
{
|
|
#if HAVE_MPI
|
|
const auto & pinfo =
|
|
boost::any_cast<const ParallelISTLInformation&>(solver_.parallelInformation());
|
|
auto operators = std::make_tuple(Opm::Reduction::makeGlobalSumFunctor<double>(),
|
|
Opm::Reduction::makeGlobalSumFunctor<double>(),
|
|
Opm::Reduction::makeGlobalSumFunctor<double>());
|
|
std::vector<double> pav_nom(p.size());
|
|
std::vector<double> pav_denom(pav_nom.size());
|
|
for (unsigned i = 0; i < p.size(); ++i) {
|
|
pav_nom[i] = p[i] * geo_.poreVolume()[i] * hydrocarbon[i];
|
|
pav_denom[i] = geo_.poreVolume()[i] * hydrocarbon[i];
|
|
}
|
|
|
|
// using ref cref to prevent copying
|
|
auto inputs = std::make_tuple(std::cref(geo_.poreVolume()),
|
|
std::cref(pav_nom), std::cref(pav_denom));
|
|
std::tuple<double, double, double> results(0.0, 0.0, 0.0);
|
|
|
|
pinfo.computeReduction(inputs, operators, results);
|
|
using std::get;
|
|
totals[5] = get<0>(results);
|
|
totals[6] = get<1>(results)/get<2>(results);
|
|
|
|
#else
|
|
// This should never happen!
|
|
OPM_THROW(std::logic_error, "HAVE_MPI should be defined if we are running in parallel");
|
|
#endif
|
|
}
|
|
|
|
return totals;
|
|
}
|
|
|
|
|
|
|
|
template <class Implementation>
|
|
void
|
|
SimulatorBase<Implementation>::outputFluidInPlace(const std::vector<double>& oip, const std::vector<double>& cip, const UnitSystem& units, const int reg)
|
|
{
|
|
std::ostringstream ss;
|
|
if (!reg) {
|
|
ss << " ===================================================\n"
|
|
<< " : Field Totals :\n";
|
|
} else {
|
|
ss << " ===================================================\n"
|
|
<< " : FIPNUM report region "
|
|
<< std::setw(2) << reg << " :\n";
|
|
}
|
|
if (units.getType() == UnitSystem::UnitType::UNIT_TYPE_METRIC) {
|
|
ss << " : PAV =" << std::setw(14) << cip[6] << " BARSA :\n"
|
|
<< std::fixed << std::setprecision(0)
|
|
<< " : PORV =" << std::setw(14) << cip[5] << " RM3 :\n";
|
|
if (!reg) {
|
|
ss << " : Pressure is weighted by hydrocarbon pore volume :\n"
|
|
<< " : Porv volumes are taken at reference conditions :\n";
|
|
}
|
|
ss << " :--------------- Oil SM3 ---------------:-- Wat SM3 --:--------------- Gas SM3 ---------------:\n";
|
|
}
|
|
if (units.getType() == UnitSystem::UnitType::UNIT_TYPE_FIELD) {
|
|
ss << " : PAV =" << std::setw(14) << cip[6] << " PSIA :\n"
|
|
<< std::fixed << std::setprecision(0)
|
|
<< " : PORV =" << std::setw(14) << cip[5] << " RB :\n";
|
|
if (!reg) {
|
|
ss << " : Pressure is weighted by hydrocarbon pore voulme :\n"
|
|
<< " : Pore volumes are taken at reference conditions :\n";
|
|
}
|
|
ss << " :--------------- Oil STB ---------------:-- Wat STB --:--------------- Gas MSCF ---------------:\n";
|
|
}
|
|
ss << " : Liquid Vapour Total : Total : Free Dissolved Total :" << "\n"
|
|
<< ":------------------------:------------------------------------------:----------------:------------------------------------------:" << "\n"
|
|
<< ":Currently in place :" << std::setw(14) << cip[1] << std::setw(14) << cip[4] << std::setw(14) << (cip[1]+cip[4]) << ":"
|
|
<< std::setw(13) << cip[0] << " :" << std::setw(14) << (cip[2]) << std::setw(14) << cip[3] << std::setw(14) << (cip[2] + cip[3]) << ":\n"
|
|
<< ":------------------------:------------------------------------------:----------------:------------------------------------------:\n"
|
|
<< ":Originally in place :" << std::setw(14) << oip[1] << std::setw(14) << oip[4] << std::setw(14) << (oip[1]+oip[4]) << ":"
|
|
<< std::setw(13) << oip[0] << " :" << std::setw(14) << oip[2] << std::setw(14) << oip[3] << std::setw(14) << (oip[2] + oip[3]) << ":\n"
|
|
<< ":========================:==========================================:================:==========================================:\n";
|
|
OpmLog::note(ss.str());
|
|
}
|
|
|
|
|
|
template <class Implementation>
|
|
void
|
|
SimulatorBase<Implementation>::
|
|
updateListEconLimited(const std::unique_ptr<Solver>& solver,
|
|
const Schedule& schedule,
|
|
const int current_step,
|
|
const Wells* wells,
|
|
const WellState& well_state,
|
|
DynamicListEconLimited& list_econ_limited) const
|
|
{
|
|
|
|
solver->model().wellModel().updateListEconLimited(schedule, current_step, wells,
|
|
well_state, list_econ_limited);
|
|
}
|
|
|
|
template <class Implementation>
|
|
void
|
|
SimulatorBase<Implementation>::
|
|
initHysteresisParams(ReservoirState& state)
|
|
{
|
|
typedef std::vector<double> VectorType;
|
|
|
|
const VectorType& somax = state.getCellData( "SOMAX" );
|
|
|
|
VectorType& pcSwMdc_ow = state.getCellData( "PCSWMDC_OW" );
|
|
VectorType& krnSwMdc_ow = state.getCellData( "KRNSWMDC_OW" );
|
|
|
|
VectorType& pcSwMdc_go = state.getCellData( "PCSWMDC_GO" );
|
|
VectorType& krnSwMdc_go = state.getCellData( "KRNSWMDC_GO" );
|
|
|
|
props_.setSatOilMax(somax);
|
|
props_.setOilWaterHystParams(pcSwMdc_ow, krnSwMdc_ow, allcells_);
|
|
props_.setGasOilHystParams(pcSwMdc_go, krnSwMdc_go, allcells_);
|
|
}
|
|
|
|
|
|
} // namespace Opm
|