2016-06-06 08:40:06 -05:00
|
|
|
/*
|
|
|
|
Copyright 2013, 2015 SINTEF ICT, Applied Mathematics.
|
|
|
|
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/>.
|
|
|
|
*/
|
|
|
|
|
|
|
|
#ifndef OPM_SIMULATORFULLYIMPLICITBLACKOIL_HEADER_INCLUDED
|
|
|
|
#define OPM_SIMULATORFULLYIMPLICITBLACKOIL_HEADER_INCLUDED
|
|
|
|
|
|
|
|
#include <opm/autodiff/SimulatorBase.hpp>
|
|
|
|
#include <opm/autodiff/NonlinearSolver.hpp>
|
|
|
|
#include <opm/autodiff/BlackoilModelEbos.hpp>
|
2016-07-12 11:47:52 -05:00
|
|
|
#include <opm/autodiff/BlackoilModelParameters.hpp>
|
2016-06-06 08:40:06 -05:00
|
|
|
|
|
|
|
namespace Opm {
|
|
|
|
|
|
|
|
class SimulatorFullyImplicitBlackoilEbos;
|
|
|
|
class StandardWells;
|
|
|
|
|
2016-07-12 11:47:52 -05:00
|
|
|
/// a simulator for the blackoil model
|
|
|
|
class SimulatorFullyImplicitBlackoilEbos
|
2016-06-06 08:40:06 -05:00
|
|
|
{
|
2016-07-12 11:47:52 -05:00
|
|
|
public:
|
2016-08-08 07:31:32 -05:00
|
|
|
typedef typename TTAG(EclFlowProblem) TypeTag;
|
2016-08-08 08:12:51 -05:00
|
|
|
typedef typename GET_PROP_TYPE(TypeTag, Simulator) Simulator;
|
|
|
|
typedef typename GET_PROP_TYPE(TypeTag, Grid) Grid;
|
2016-06-06 08:40:06 -05:00
|
|
|
typedef WellStateFullyImplicitBlackoil WellState;
|
|
|
|
typedef BlackoilState ReservoirState;
|
|
|
|
typedef BlackoilOutputWriter OutputWriter;
|
2016-08-08 07:41:52 -05:00
|
|
|
typedef BlackoilModelEbos Model;
|
2016-07-12 11:47:52 -05:00
|
|
|
typedef BlackoilModelParameters ModelParameters;
|
2016-06-06 08:40:06 -05:00
|
|
|
typedef NonlinearSolver<Model> Solver;
|
|
|
|
typedef StandardWells WellModel;
|
|
|
|
|
2016-07-12 11:47:52 -05:00
|
|
|
|
|
|
|
/// Initialise from parameters and objects to observe.
|
|
|
|
/// \param[in] param parameters, this class accepts the following:
|
|
|
|
/// parameter (default) effect
|
|
|
|
/// -----------------------------------------------------------
|
|
|
|
/// output (true) write output to files?
|
|
|
|
/// output_dir ("output") output directoty
|
|
|
|
/// output_interval (1) output every nth step
|
|
|
|
/// nl_pressure_residual_tolerance (0.0) pressure solver residual tolerance (in Pascal)
|
|
|
|
/// nl_pressure_change_tolerance (1.0) pressure solver change tolerance (in Pascal)
|
|
|
|
/// nl_pressure_maxiter (10) max nonlinear iterations in pressure
|
|
|
|
/// nl_maxiter (30) max nonlinear iterations in transport
|
|
|
|
/// nl_tolerance (1e-9) transport solver absolute residual tolerance
|
|
|
|
/// num_transport_substeps (1) number of transport steps per pressure step
|
|
|
|
/// use_segregation_split (false) solve for gravity segregation (if false,
|
|
|
|
/// segregation is ignored).
|
|
|
|
///
|
|
|
|
/// \param[in] geo derived geological properties
|
|
|
|
/// \param[in] props fluid and rock properties
|
|
|
|
/// \param[in] rock_comp_props if non-null, rock compressibility properties
|
|
|
|
/// \param[in] linsolver linear solver
|
|
|
|
/// \param[in] gravity if non-null, gravity vector
|
|
|
|
/// \param[in] has_disgas true for dissolved gas option
|
|
|
|
/// \param[in] has_vapoil true for vaporized oil option
|
|
|
|
/// \param[in] eclipse_state the object which represents an internalized ECL deck
|
|
|
|
/// \param[in] output_writer
|
|
|
|
/// \param[in] threshold_pressures_by_face if nonempty, threshold pressures that inhibit flow
|
2016-08-08 08:26:09 -05:00
|
|
|
SimulatorFullyImplicitBlackoilEbos(Simulator& ebosSimulator,
|
|
|
|
const parameter::ParameterGroup& param,
|
2016-07-12 11:47:52 -05:00
|
|
|
DerivedGeology& geo,
|
|
|
|
BlackoilPropsAdInterface& 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,
|
|
|
|
BlackoilOutputWriter& output_writer,
|
|
|
|
const std::vector<double>& threshold_pressures_by_face)
|
2016-08-08 08:26:09 -05:00
|
|
|
: ebosSimulator_(ebosSimulator),
|
|
|
|
param_(param),
|
2016-07-12 11:47:52 -05:00
|
|
|
model_param_(param),
|
|
|
|
solver_param_(param),
|
|
|
|
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)),
|
|
|
|
output_writer_(output_writer),
|
|
|
|
threshold_pressures_by_face_(threshold_pressures_by_face),
|
|
|
|
is_parallel_run_( false )
|
|
|
|
{
|
2016-08-08 07:31:32 -05:00
|
|
|
|
2016-08-08 08:12:51 -05:00
|
|
|
// Misc init.
|
|
|
|
const int num_cells = AutoDiffGrid::numCells(grid());
|
|
|
|
allcells_.resize(num_cells);
|
|
|
|
for (int cell = 0; cell < num_cells; ++cell) {
|
|
|
|
allcells_[cell] = cell;
|
|
|
|
}
|
|
|
|
|
|
|
|
rateConverter_.reset(new RateConverterType(props_, std::vector<int>(AutoDiffGrid::numCells(grid()), 0)));
|
|
|
|
|
2016-07-12 11:47:52 -05:00
|
|
|
#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
|
|
|
|
}
|
|
|
|
|
|
|
|
/// Run the simulation.
|
|
|
|
/// This will run succesive timesteps until timer.done() is true. It will
|
|
|
|
/// modify the reservoir and well states.
|
|
|
|
/// \param[in,out] timer governs the requested reporting timesteps
|
|
|
|
/// \param[in,out] state state of reservoir: pressure, fluxes
|
|
|
|
/// \return simulation report, with timing data
|
|
|
|
SimulatorReport run(SimulatorTimer& timer,
|
|
|
|
ReservoirState& state)
|
|
|
|
{
|
|
|
|
WellState prev_well_state;
|
|
|
|
|
|
|
|
if (output_writer_.isRestart()) {
|
|
|
|
// This is a restart, populate WellState and ReservoirState state objects from restart file
|
2016-08-08 08:12:51 -05:00
|
|
|
output_writer_.initFromRestartFile(props_.phaseUsage(), props_.permeability(), grid(), state, prev_well_state);
|
|
|
|
initHydroCarbonState(state, props_.phaseUsage(), Opm::UgGridHelpers::numCells(grid()), has_disgas_, has_vapoil_);
|
2016-07-12 11:47:52 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
// Create timers and file for writing timing info.
|
|
|
|
Opm::time::StopWatch solver_timer;
|
|
|
|
double stime = 0.0;
|
|
|
|
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(tstep_filename.c_str());
|
|
|
|
|
2016-08-08 08:26:09 -05:00
|
|
|
const auto& schedule = eclState()->getSchedule();
|
2016-07-12 11:47:52 -05:00
|
|
|
const auto& events = schedule->getEvents();
|
|
|
|
|
|
|
|
// adaptive time stepping
|
|
|
|
std::unique_ptr< AdaptiveTimeStepping > adaptiveTimeStepping;
|
|
|
|
if( param_.getDefault("timestep.adaptive", true ) )
|
|
|
|
{
|
|
|
|
adaptiveTimeStepping.reset( new AdaptiveTimeStepping( param_, terminal_output_ ) );
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
2016-08-08 08:12:51 -05:00
|
|
|
output_writer_.writeInit( geo_.simProps(grid()) , geo_.nonCartesianConnections( ) );
|
2016-07-12 11:47:52 -05:00
|
|
|
|
|
|
|
std::string restorefilename = param_.getDefault("restorefile", std::string("") );
|
|
|
|
if( ! restorefilename.empty() )
|
|
|
|
{
|
|
|
|
// -1 means that we'll take the last report step that was written
|
|
|
|
const int desiredRestoreStep = param_.getDefault("restorestep", int(-1) );
|
|
|
|
|
|
|
|
output_writer_.restore( timer,
|
|
|
|
state,
|
|
|
|
prev_well_state,
|
|
|
|
restorefilename,
|
|
|
|
desiredRestoreStep );
|
|
|
|
}
|
|
|
|
|
|
|
|
unsigned int totalLinearizations = 0;
|
|
|
|
unsigned int totalNonlinearIterations = 0;
|
|
|
|
unsigned int totalLinearIterations = 0;
|
|
|
|
bool is_well_potentials_computed = param_.getDefault("compute_well_potentials", false );
|
|
|
|
std::vector<double> well_potentials;
|
|
|
|
DynamicListEconLimited dynamic_list_econ_limited;
|
|
|
|
|
|
|
|
// 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.
|
2016-08-08 08:26:09 -05:00
|
|
|
WellsManager wells_manager(eclState(),
|
2016-07-12 11:47:52 -05:00
|
|
|
timer.currentStepNum(),
|
2016-08-08 08:12:51 -05:00
|
|
|
Opm::UgGridHelpers::numCells(grid()),
|
|
|
|
Opm::UgGridHelpers::globalCell(grid()),
|
|
|
|
Opm::UgGridHelpers::cartDims(grid()),
|
|
|
|
Opm::UgGridHelpers::dimensions(grid()),
|
|
|
|
Opm::UgGridHelpers::cell2Faces(grid()),
|
|
|
|
Opm::UgGridHelpers::beginFaceCentroids(grid()),
|
2016-07-12 11:47:52 -05:00
|
|
|
props_.permeability(),
|
|
|
|
dynamic_list_econ_limited,
|
|
|
|
is_parallel_run_,
|
|
|
|
well_potentials);
|
|
|
|
const Wells* wells = wells_manager.c_wells();
|
|
|
|
WellState well_state;
|
|
|
|
well_state.init(wells, state, prev_well_state);
|
|
|
|
|
|
|
|
// give the polymer and surfactant simulators the chance to do their stuff
|
|
|
|
handleAdditionalWellInflow(timer, wells_manager, well_state, wells);
|
|
|
|
|
|
|
|
// write the inital state at the report stage
|
|
|
|
if (timer.initialStep()) {
|
|
|
|
output_writer_.writeTimeStep( timer, state, well_state );
|
|
|
|
}
|
|
|
|
|
|
|
|
// Compute reservoir volumes for RESV controls.
|
|
|
|
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);
|
|
|
|
|
|
|
|
auto solver = createSolver(well_model);
|
|
|
|
|
|
|
|
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()
|
|
|
|
<< "\n";
|
|
|
|
OpmLog::info(step_msg.str());
|
|
|
|
}
|
|
|
|
|
2016-08-08 07:58:25 -05:00
|
|
|
solver->model().beginReportStep();
|
|
|
|
|
2016-07-12 11:47:52 -05:00
|
|
|
// 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 ) {
|
|
|
|
adaptiveTimeStepping->step( timer, *solver, state, well_state, output_writer_ );
|
|
|
|
}
|
|
|
|
else {
|
|
|
|
// solve for complete report step
|
|
|
|
solver->step(timer, state, well_state);
|
2016-08-14 14:44:47 -05:00
|
|
|
|
|
|
|
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());
|
|
|
|
}
|
2016-07-12 11:47:52 -05:00
|
|
|
}
|
|
|
|
|
2016-08-08 07:58:25 -05:00
|
|
|
solver->model().endReportStep();
|
|
|
|
|
2016-07-12 11:47:52 -05:00
|
|
|
// 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)
|
|
|
|
DeckConstPtr miniDeck = schedule->getModifierDeck(nextTimeStepIdx);
|
|
|
|
eclState()->applyModifierDeck(*miniDeck);
|
2016-08-08 08:26:09 -05:00
|
|
|
geo_.update(grid(), props_, eclState(), gravity_);
|
2016-07-12 11:47:52 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
// take time that was used to solve system for this reportStep
|
|
|
|
solver_timer.stop();
|
|
|
|
|
|
|
|
// accumulate the number of nonlinear and linear Iterations
|
|
|
|
totalLinearizations += solver->linearizations();
|
|
|
|
totalNonlinearIterations += solver->nonlinearIterations();
|
|
|
|
totalLinearIterations += solver->linearIterations();
|
|
|
|
|
|
|
|
// Report timing.
|
|
|
|
const double st = solver_timer.secsSinceStart();
|
|
|
|
|
|
|
|
// accumulate total time
|
|
|
|
stime += st;
|
|
|
|
|
|
|
|
if ( terminal_output_ )
|
|
|
|
{
|
|
|
|
std::string msg;
|
|
|
|
msg = "Fully implicit solver took: " + std::to_string(st) + " seconds. Total solver time taken: " + std::to_string(stime) + " seconds.";
|
|
|
|
OpmLog::note(msg);
|
|
|
|
}
|
|
|
|
|
|
|
|
if ( output_writer_.output() ) {
|
|
|
|
SimulatorReport step_report;
|
|
|
|
step_report.pressure_time = st;
|
|
|
|
step_report.total_time = step_timer.secsSinceStart();
|
|
|
|
step_report.reportParam(tstep_os);
|
|
|
|
}
|
|
|
|
|
|
|
|
// Increment timer, remember well state.
|
|
|
|
++timer;
|
|
|
|
|
|
|
|
// write simulation state at the report stage
|
|
|
|
output_writer_.writeTimeStep( timer, state, well_state );
|
|
|
|
|
|
|
|
prev_well_state = well_state;
|
|
|
|
// The well potentials are only computed if they are needed
|
|
|
|
// For now thay are only used to determine default guide rates for group controlled wells
|
|
|
|
if ( is_well_potentials_computed ) {
|
|
|
|
computeWellPotentials(wells, well_state, well_potentials);
|
|
|
|
}
|
|
|
|
|
2016-08-08 08:26:09 -05:00
|
|
|
updateListEconLimited(solver, eclState()->getSchedule(), timer.currentStepNum(), wells,
|
2016-07-12 11:47:52 -05:00
|
|
|
well_state, dynamic_list_econ_limited);
|
|
|
|
}
|
|
|
|
|
|
|
|
// Stop timer and create timing report
|
|
|
|
total_timer.stop();
|
|
|
|
SimulatorReport report;
|
|
|
|
report.pressure_time = stime;
|
|
|
|
report.transport_time = 0.0;
|
|
|
|
report.total_time = total_timer.secsSinceStart();
|
|
|
|
report.total_linearizations = totalLinearizations;
|
|
|
|
report.total_newton_iterations = totalNonlinearIterations;
|
|
|
|
report.total_linear_iterations = totalLinearIterations;
|
|
|
|
return report;
|
|
|
|
}
|
|
|
|
|
2016-08-08 08:12:51 -05:00
|
|
|
const Grid& grid() const
|
2016-08-08 08:26:09 -05:00
|
|
|
{ return ebosSimulator_.gridManager().grid(); }
|
2016-08-08 08:12:51 -05:00
|
|
|
|
2016-07-12 11:47:52 -05:00
|
|
|
protected:
|
|
|
|
void handleAdditionalWellInflow(SimulatorTimer& timer,
|
|
|
|
WellsManager& wells_manager,
|
|
|
|
WellState& well_state,
|
|
|
|
const Wells* wells)
|
|
|
|
{ }
|
|
|
|
|
|
|
|
std::unique_ptr<Solver> createSolver(const WellModel& well_model)
|
|
|
|
{
|
2016-08-08 08:26:09 -05:00
|
|
|
auto model = std::unique_ptr<Model>(new Model(ebosSimulator_,
|
2016-08-08 07:31:32 -05:00
|
|
|
model_param_,
|
2016-07-12 11:47:52 -05:00
|
|
|
props_,
|
|
|
|
geo_,
|
|
|
|
rock_comp_props_,
|
|
|
|
well_model,
|
|
|
|
solver_,
|
|
|
|
terminal_output_));
|
|
|
|
|
|
|
|
return std::unique_ptr<Solver>(new Solver(solver_param_, std::move(model)));
|
|
|
|
}
|
|
|
|
|
|
|
|
void computeRESV(const std::size_t step,
|
|
|
|
const Wells* wells,
|
|
|
|
const BlackoilState& x,
|
|
|
|
WellState& xw)
|
|
|
|
{
|
|
|
|
typedef SimFIBODetails::WellMap WellMap;
|
|
|
|
|
2016-08-08 08:26:09 -05:00
|
|
|
const auto w_ecl = eclState()->getSchedule()->getWells(step);
|
2016-07-12 11:47:52 -05:00
|
|
|
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.
|
2016-08-08 08:12:51 -05:00
|
|
|
rateConverter_->defineState(x, boost::any_cast<const ParallelISTLInformation&>(solver_.parallelInformation()));
|
2016-07-12 11:47:52 -05:00
|
|
|
}
|
|
|
|
}
|
|
|
|
else
|
|
|
|
#endif
|
|
|
|
{
|
|
|
|
if ( global_number_resv_wells )
|
|
|
|
{
|
2016-08-08 08:12:51 -05:00
|
|
|
rateConverter_->defineState(x);
|
2016-07-12 11:47:52 -05:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
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.
|
2016-08-08 08:12:51 -05:00
|
|
|
rateConverter_->calcCoeff(prates, fipreg, distr);
|
2016-07-12 11:47:52 -05:00
|
|
|
|
|
|
|
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.
|
2016-08-08 08:12:51 -05:00
|
|
|
rateConverter_->calcCoeff(hrates, fipreg, distr);
|
2016-07-12 11:47:52 -05:00
|
|
|
|
|
|
|
// 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.");
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
void computeWellPotentials(const Wells* wells,
|
|
|
|
const WellState& xw,
|
|
|
|
std::vector<double>& well_potentials)
|
|
|
|
{
|
|
|
|
const int nw = wells->number_of_wells;
|
|
|
|
const int np = wells->number_of_phases;
|
|
|
|
well_potentials.clear();
|
|
|
|
well_potentials.resize(nw*np,0.0);
|
|
|
|
for (int w = 0; w < nw; ++w) {
|
|
|
|
for (int perf = wells->well_connpos[w]; perf < wells->well_connpos[w + 1]; ++perf) {
|
|
|
|
for (int phase = 0; phase < np; ++phase) {
|
|
|
|
well_potentials[w*np + phase] += xw.wellPotentials()[perf*np + phase];
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
void updateListEconLimited(const std::unique_ptr<Solver>& solver,
|
|
|
|
ScheduleConstPtr 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);
|
|
|
|
}
|
|
|
|
|
2016-08-08 08:26:09 -05:00
|
|
|
EclipseStateConstPtr eclState() const
|
|
|
|
{ return ebosSimulator_.gridManager().eclState(); }
|
2016-07-12 11:47:52 -05:00
|
|
|
|
2016-08-08 08:26:09 -05:00
|
|
|
EclipseStatePtr eclState()
|
|
|
|
{ return ebosSimulator_.gridManager().eclState(); }
|
2016-07-12 11:47:52 -05:00
|
|
|
|
|
|
|
// Data.
|
2016-08-08 08:26:09 -05:00
|
|
|
Simulator& ebosSimulator_;
|
2016-08-08 07:31:32 -05:00
|
|
|
|
2016-07-12 11:47:52 -05:00
|
|
|
typedef RateConverter::
|
|
|
|
SurfaceToReservoirVoidage< BlackoilPropsAdInterface,
|
|
|
|
std::vector<int> > RateConverterType;
|
|
|
|
typedef typename Solver::SolverParameters SolverParameters;
|
|
|
|
|
|
|
|
const parameter::ParameterGroup param_;
|
|
|
|
ModelParameters model_param_;
|
|
|
|
SolverParameters solver_param_;
|
|
|
|
|
|
|
|
// Observed objects.
|
|
|
|
BlackoilPropsAdInterface& props_;
|
|
|
|
const RockCompressibility* rock_comp_props_;
|
|
|
|
const double* gravity_;
|
|
|
|
// Solvers
|
|
|
|
DerivedGeology& geo_;
|
|
|
|
NewtonIterationBlackoilInterface& solver_;
|
|
|
|
// Misc. data
|
|
|
|
std::vector<int> allcells_;
|
|
|
|
const bool has_disgas_;
|
|
|
|
const bool has_vapoil_;
|
|
|
|
bool terminal_output_;
|
|
|
|
// output_writer
|
|
|
|
OutputWriter& output_writer_;
|
2016-08-08 08:12:51 -05:00
|
|
|
std::unique_ptr<RateConverterType> rateConverter_;
|
2016-07-12 11:47:52 -05:00
|
|
|
// Threshold pressures.
|
|
|
|
std::vector<double> threshold_pressures_by_face_;
|
|
|
|
// Whether this a parallel simulation or not
|
|
|
|
bool is_parallel_run_;
|
|
|
|
|
2016-06-06 08:40:06 -05:00
|
|
|
};
|
|
|
|
|
|
|
|
} // namespace Opm
|
|
|
|
|
|
|
|
#endif // OPM_SIMULATORFULLYIMPLICITBLACKOIL_HEADER_INCLUDED
|