making rate_converter to be reference to the one in Simulator

keeping the const property in the Well Model.
This commit is contained in:
Kai Bao 2017-08-10 11:20:09 +02:00
parent 1950052684
commit c59aa9127e
4 changed files with 29 additions and 40 deletions

View File

@ -155,6 +155,7 @@ namespace Opm {
BlackoilModelEbos(Simulator& ebosSimulator,
const ModelParameters& param,
const StandardWellsDense<TypeTag>& well_model,
RateConverterType& rate_converter,
const NewtonIterationBlackoilInterface& linsolver,
const bool terminal_output
)
@ -173,7 +174,7 @@ namespace Opm {
, param_( param )
, well_model_ (well_model)
, terminal_output_ (terminal_output)
, rate_converter_(wellModel().rateConverter())
, rate_converter_(rate_converter)
, current_relaxation_(1.0)
, dx_old_(AutoDiffGrid::numCells(grid_))
, isBeginReportStep_(false)
@ -1505,7 +1506,7 @@ namespace Opm {
long int global_nc_;
// rate converter between the surface volume rates and reservoir voidage rates
RateConverterType* rate_converter_;
RateConverterType& rate_converter_;
std::vector<std::vector<double>> residual_norms_history_;
double current_relaxation_;
@ -1660,7 +1661,7 @@ namespace Opm {
global_number_wells = info.communicator().sum(global_number_wells);
if ( global_number_wells )
{
rate_converter_->defineState(reservoir_state, boost::any_cast<const ParallelISTLInformation&>(istlSolver_->parallelInformation()));
rate_converter_.defineState(reservoir_state, boost::any_cast<const ParallelISTLInformation&>(istlSolver_->parallelInformation()));
}
}
else
@ -1668,7 +1669,7 @@ namespace Opm {
{
if ( global_number_wells )
{
rate_converter_->defineState(reservoir_state);
rate_converter_.defineState(reservoir_state);
}
}
}

View File

@ -65,6 +65,7 @@ public:
typedef BlackoilModelParameters ModelParameters;
typedef NonlinearSolver<Model> Solver;
typedef StandardWellsDense<TypeTag> WellModel;
typedef RateConverter::SurfaceToReservoirVoidage<FluidSystem, std::vector<int> > RateConverterType;
/// Initialise from parameters and objects to observe.
@ -108,14 +109,10 @@ public:
has_vapoil_(has_vapoil),
terminal_output_(param.getDefault("output_terminal", true)),
output_writer_(output_writer),
rateConverter_(createRateConverter_()),
defunct_well_names_( defunct_well_names ),
is_parallel_run_( false )
{
extractLegacyCellPvtRegionIndex_();
rateConverter_.reset(new RateConverterType(phaseUsage_,
legacyCellPvtRegionIdx_.data(),
AutoDiffGrid::numCells(grid()),
std::vector<int>(AutoDiffGrid::numCells(grid()), 0)));
#if HAVE_MPI
if ( solver_.parallelInformation().type() == typeid(ParallelISTLInformation) )
@ -254,8 +251,8 @@ public:
solver_timer.start();
const auto& wells_ecl = eclState().getSchedule().getWells(timer.currentStepNum());
WellModel well_model(wells, &(wells_manager.wellCollection()), wells_ecl, model_param_, terminal_output_,
timer.currentStepNum());
WellModel well_model(wells, &(wells_manager.wellCollection()), wells_ecl, model_param_,
rateConverter_, terminal_output_, timer.currentStepNum());
auto solver = createSolver(well_model);
@ -438,12 +435,12 @@ protected:
gravity,
legacyDepth_,
legacyPoreVolume_,
rateConverter_.get(),
globalNumCells,
grid());
auto model = std::unique_ptr<Model>(new Model(ebosSimulator_,
model_param_,
well_model,
rateConverter_,
solver_,
terminal_output_));
@ -476,7 +473,7 @@ protected:
// 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()));
rateConverter_.defineState(x, boost::any_cast<const ParallelISTLInformation&>(solver_.parallelInformation()));
}
}
else
@ -484,7 +481,7 @@ protected:
{
if ( global_number_resv_wells )
{
rateConverter_->defineState(x);
rateConverter_.defineState(x);
}
}
@ -523,7 +520,7 @@ protected:
}
const int fipreg = 0; // Hack. Ignore FIP regions.
rateConverter_->calcCoeff(prates, fipreg, distr);
rateConverter_.calcCoeff(prates, fipreg, distr);
well_controls_iset_distr(ctrl, rctrl, & distr[0]);
}
@ -545,7 +542,7 @@ protected:
SimFIBODetails::historyRates(pu, p, hrates);
const int fipreg = 0; // Hack. Ignore FIP regions.
rateConverter_->calcCoeff(hrates, fipreg, distr);
rateConverter_.calcCoeff(hrates, fipreg, distr);
// WCONHIST/RESV target is sum of all
// observed phase rates translated to
@ -853,8 +850,14 @@ protected:
}
}
RateConverterType createRateConverter_() {
extractLegacyCellPvtRegionIndex_();
RateConverterType rate_converter(phaseUsage_,
legacyCellPvtRegionIdx_.data(),
AutoDiffGrid::numCells(grid()),
std::vector<int>(AutoDiffGrid::numCells(grid()), 0));
return rate_converter;
}
// Data.
@ -863,7 +866,6 @@ protected:
std::vector<int> legacyCellPvtRegionIdx_;
std::vector<double> legacyPoreVolume_;
std::vector<double> legacyDepth_;
typedef RateConverter::SurfaceToReservoirVoidage<FluidSystem, std::vector<int> > RateConverterType;
typedef typename Solver::SolverParameters SolverParameters;
SimulatorReport failureReport_;
@ -881,7 +883,7 @@ protected:
bool terminal_output_;
// output_writer
OutputWriter& output_writer_;
std::unique_ptr<RateConverterType> rateConverter_;
RateConverterType rateConverter_;
// The names of wells that should be defunct
// (e.g. in a parallel run when they are handeled by
// a different process)

View File

@ -113,6 +113,7 @@ enum WellVariablePositions {
WellCollection* well_collection,
const std::vector< const Well* >& wells_ecl,
const ModelParameters& param,
const RateConverterType& rate_converter,
const bool terminal_output,
const int current_index);
@ -121,7 +122,6 @@ enum WellVariablePositions {
const double gravity_arg,
const std::vector<double>& depth_arg,
const std::vector<double>& pv_arg,
RateConverterType* rate_converter,
long int global_nc,
const Grid& grid);
@ -291,8 +291,6 @@ enum WellVariablePositions {
void applyVREPGroupControl(WellState& well_state) const;
RateConverterType* rateConverter() const;
protected:
bool wells_active_;
const Wells* wells_;
@ -311,7 +309,7 @@ enum WellVariablePositions {
std::vector<bool> active_;
const VFPProperties* vfp_properties_;
double gravity_;
RateConverterType* rate_converter_;
const RateConverterType& rate_converter_;
// The efficiency factor for each connection. It is specified based on wells and groups,
// We calculate the factor for each connection for the computation of contributions to the mass balance equations.

View File

@ -9,6 +9,7 @@ namespace Opm {
WellCollection* well_collection,
const std::vector< const Well* >& wells_ecl,
const ModelParameters& param,
const RateConverterType& rate_converter,
const bool terminal_output,
const int current_timeIdx)
: wells_active_(wells_arg!=nullptr)
@ -20,6 +21,7 @@ namespace Opm {
, has_solvent_(GET_PROP_VALUE(TypeTag, EnableSolvent))
, has_polymer_(GET_PROP_VALUE(TypeTag, EnablePolymer))
, current_timeIdx_(current_timeIdx)
, rate_converter_(rate_converter)
, well_perforation_efficiency_factors_((wells_!=nullptr ? wells_->well_connpos[wells_->number_of_wells] : 0), 1.0)
, well_perforation_densities_( wells_ ? wells_arg->well_connpos[wells_arg->number_of_wells] : 0)
, well_perforation_pressure_diffs_( wells_ ? wells_arg->well_connpos[wells_arg->number_of_wells] : 0)
@ -46,7 +48,6 @@ namespace Opm {
const double gravity_arg,
const std::vector<double>& depth_arg,
const std::vector<double>& pv_arg,
RateConverterType* rate_converter,
long int global_nc,
const Grid& grid)
{
@ -62,7 +63,6 @@ namespace Opm {
gravity_ = gravity_arg;
cell_depths_ = extractPerfData(depth_arg);
pv_ = pv_arg;
rate_converter_ = rate_converter;
calculateEfficiencyFactors();
@ -144,18 +144,6 @@ namespace Opm {
template<typename TypeTag>
typename StandardWellsDense<TypeTag>::RateConverterType*
StandardWellsDense<TypeTag>::
rateConverter() const
{
return rate_converter_;
}
template<typename TypeTag>
SimulatorReport
StandardWellsDense<TypeTag>::
@ -2123,7 +2111,7 @@ namespace Opm {
// the average hydrocarbon conditions of the whole field will be used
const int fipreg = 0; // Not considering FIP for the moment.
rate_converter_->calcCoeff(well_rates, fipreg, convert_coeff);
rate_converter_.calcCoeff(well_rates, fipreg, convert_coeff);
well_voidage_rates[w] = std::inner_product(well_rates.begin(), well_rates.end(),
convert_coeff.begin(), 0.0);
} else {
@ -2134,7 +2122,7 @@ namespace Opm {
well_rates.begin());
// the average hydrocarbon conditions of the whole field will be used
const int fipreg = 0; // Not considering FIP for the moment.
rate_converter_->calcCoeff(well_rates, fipreg, convert_coeff);
rate_converter_.calcCoeff(well_rates, fipreg, convert_coeff);
std::copy(convert_coeff.begin(), convert_coeff.end(),
voidage_conversion_coeffs.begin() + np * w);
}